Ignore:
Timestamp:
Feb 2, 2012, 1:58:12 PM (14 years ago)
Author:
Julian Iseringhausen <isering@…>
Children:
32ff22
Parents:
01be70
Message:

Parallel performance update.

git-svn-id: https://svn.version.fz-juelich.de/scafacos/trunk@1314 5161e1c8-67bf-11de-9fd5-51895aff932f

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/smoother/gsrb.cpp

    r01be70 r894a5f  
    1212#endif
    1313
     14#include <sstream>
     15
    1416#include "base/discretization.hpp"
    15 #include "base/timer.hpp"
    1617#include "base/stencil.hpp"
    1718#include "comm/comm.hpp"
     
    2021
    2122using namespace VMG;
     23
     24static inline void ComputePartial(Grid& sol, Grid& rhs, const Stencil& mat,
     25                                  const Index& begin, const Index& end,
     26                                  const vmg_float& prefactor, const vmg_float& diag_inv,
     27                                  const int& off)
     28{
     29  int i,j,k;
     30  vmg_float temp;
     31  Stencil::iterator iter;
     32
     33  for (i=begin.X(); i<end.X(); ++i)
     34    for (j=begin.Y(); j<end.Y(); ++j)
     35      for (k=begin.Z()+((i+j+off)%2); k<end.Z(); k+=2) {
     36
     37        temp = prefactor * rhs.GetVal(i,j,k);
     38
     39        for (iter=mat.begin(); iter!=mat.end(); ++iter)
     40          temp -= iter->Val() * sol.GetVal(i+iter->Disp().X(),
     41                                           j+iter->Disp().Y(),
     42                                           k+iter->Disp().Z());
     43
     44        sol(i,j,k) = temp * diag_inv;
     45
     46      }
     47}
    2248
    2349void GaussSeidelRB::Compute(Grid& sol, Grid& rhs)
     
    3258#endif
    3359
    34   Timer::Start("SmootherWithoutCommunication");
     60  const Stencil& mat = MG::GetDiscretization()->GetStencil();
     61  const vmg_float prefactor_inv = 1.0 / MG::GetDiscretization()->OperatorPrefactor(sol);
     62  const vmg_float diag_inv = 1.0 / mat.GetDiag();
     63  const int off = rhs.Global().BeginLocal().Sum() - rhs.Local().HasHalo1().Sum();
     64  Comm& comm = *MG::GetComm();
    3565
    36   Index i;
    37   vmg_float temp;
    38   const Stencil& A = MG::GetDiscretization()->GetStencil();
    39   const vmg_float prefactor_inv = 1.0 / MG::GetDiscretization()->OperatorPrefactor(sol);
    40   const vmg_float diag_inv = 1.0 / A.GetDiag();
    41   const int off = rhs.Global().BeginLocal().Sum() - rhs.Local().HasHalo1().Sum();
     66  /*
     67   * Compute first halfstep
     68   */
    4269
    43   for (i.X()=rhs.Local().Begin().X(); i.X()<rhs.Local().End().X(); ++i.X())
    44     for (i.Y()=rhs.Local().Begin().Y(); i.Y()<rhs.Local().End().Y(); ++i.Y())
    45       for (i.Z()=rhs.Local().Begin().Z() + ((i.X()+i.Y()+off+1) % 2); i.Z()<rhs.Local().End().Z(); i.Z()+=2) {
     70  // Start asynchronous communication
     71  comm.CommToGhostsAsyncStart(sol);
    4672
    47         temp = prefactor_inv * rhs.GetVal(i);
     73  // Smooth part not depending on ghost cells
     74  ComputePartial(sol, rhs, mat,
     75                 rhs.Local().Begin()+1, rhs.Local().End()-1,
     76                 prefactor_inv, diag_inv, off+1);
    4877
    49         for (Stencil::iterator iter=A.begin(); iter!=A.end(); ++iter)
    50           temp -= iter->Val() * sol.GetVal(i + iter->Disp());
     78  // Finish asynchronous communication
     79  comm.CommToGhostsAsyncFinish();
    5180
    52         sol(i) = temp * diag_inv;
     81  /*
     82   * Smooth near boundary cells
     83   */
    5384
    54       }
     85  ComputePartial(sol, rhs, mat,
     86                 rhs.Local().Begin(),
     87                 Index(rhs.Local().Begin().X()+1, rhs.Local().End().Y(), rhs.Local().End().Z()),
     88                 prefactor_inv, diag_inv, off+1);
    5589
    56   Timer::Stop("SmootherWithoutCommunication");
     90  ComputePartial(sol, rhs, mat,
     91                 Index(rhs.Local().End().X()-1, rhs.Local().Begin().Y(), rhs.Local().Begin().Z()),
     92                 rhs.Local().End(),
     93                 prefactor_inv, diag_inv, off+1);
    5794
    58   MG::GetComm()->CommToGhosts(sol);
     95  ComputePartial(sol, rhs, mat,
     96                 Index(rhs.Local().Begin().X()+1, rhs.Local().Begin().Y(), rhs.Local().Begin().Z()),
     97                 Index(rhs.Local().End().X()-1, rhs.Local().Begin().Y()+1, rhs.Local().End().Z()),
     98                 prefactor_inv, diag_inv, off+1);
    5999
    60   Timer::Start("SmootherWithoutCommunication");
     100  ComputePartial(sol, rhs, mat,
     101                 Index(rhs.Local().Begin().X()+1, rhs.Local().End().Y()-1, rhs.Local().Begin().Z()),
     102                 Index(rhs.Local().End().X()-1, rhs.Local().End().Y(), rhs.Local().End().Z()),
     103                 prefactor_inv, diag_inv, off+1);
    61104
    62   for (i.X()=rhs.Local().Begin().X(); i.X()<rhs.Local().End().X(); ++i.X())
    63     for (i.Y()=rhs.Local().Begin().Y(); i.Y()<rhs.Local().End().Y(); ++i.Y())
    64       for (i.Z()=rhs.Local().Begin().Z() + ((i.X()+i.Y()+off) % 2); i.Z()<rhs.Local().End().Z(); i.Z()+=2) {
     105  ComputePartial(sol, rhs, mat,
     106                 Index(rhs.Local().Begin().X()+1, rhs.Local().Begin().Y()+1, rhs.Local().Begin().Z()),
     107                 Index(rhs.Local().End().X()-1, rhs.Local().End().Y()-1, rhs.Local().Begin().Z()+1),
     108                 prefactor_inv, diag_inv, off+1);
    65109
    66         temp = prefactor_inv * rhs.GetVal(i);
     110  ComputePartial(sol, rhs, mat,
     111                 Index(rhs.Local().Begin().X()+1, rhs.Local().Begin().Y()+1, rhs.Local().End().Z()-1),
     112                 Index(rhs.Local().End().X()-1, rhs.Local().End().Y()-1, rhs.Local().End().Z()),
     113                 prefactor_inv, diag_inv, off+1);
    67114
    68         for (Stencil::iterator iter=A.begin(); iter!=A.end(); ++iter)
    69           temp -= iter->Val() * sol.GetVal(i + iter->Disp());
     115  /*
     116   * Compute second halfstep
     117   */
    70118
    71         sol(i) = temp * diag_inv;
     119  // Start asynchronous communication
     120  comm.CommToGhostsAsyncStart(sol);
    72121
    73       }
     122  // Smooth part not depending on ghost cells
     123  ComputePartial(sol, rhs, mat,
     124                 rhs.Local().Begin()+1, rhs.Local().End()-1,
     125                 prefactor_inv, diag_inv, off);
    74126
    75   Timer::Stop("SmootherWithoutCommunication");
     127  // Finish asynchronous communication
     128  comm.CommToGhostsAsyncFinish();
     129
     130  /*
     131   * Smooth near boundary cells
     132   */
     133
     134  ComputePartial(sol, rhs, mat,
     135                 rhs.Local().Begin(),
     136                 Index(rhs.Local().Begin().X()+1, rhs.Local().End().Y(), rhs.Local().End().Z()),
     137                 prefactor_inv, diag_inv, off);
     138
     139  ComputePartial(sol, rhs, mat,
     140                 Index(rhs.Local().End().X()-1, rhs.Local().Begin().Y(), rhs.Local().Begin().Z()),
     141                 rhs.Local().End(),
     142                 prefactor_inv, diag_inv, off);
     143
     144  ComputePartial(sol, rhs, mat,
     145                 Index(rhs.Local().Begin().X()+1, rhs.Local().Begin().Y(), rhs.Local().Begin().Z()),
     146                 Index(rhs.Local().End().X()-1, rhs.Local().Begin().Y()+1, rhs.Local().End().Z()),
     147                 prefactor_inv, diag_inv, off);
     148
     149  ComputePartial(sol, rhs, mat,
     150                 Index(rhs.Local().Begin().X()+1, rhs.Local().End().Y()-1, rhs.Local().Begin().Z()),
     151                 Index(rhs.Local().End().X()-1, rhs.Local().End().Y(), rhs.Local().End().Z()),
     152                 prefactor_inv, diag_inv, off);
     153
     154  ComputePartial(sol, rhs, mat,
     155                 Index(rhs.Local().Begin().X()+1, rhs.Local().Begin().Y()+1, rhs.Local().Begin().Z()),
     156                 Index(rhs.Local().End().X()-1, rhs.Local().End().Y()-1, rhs.Local().Begin().Z()+1),
     157                 prefactor_inv, diag_inv, off);
     158
     159  ComputePartial(sol, rhs, mat,
     160                 Index(rhs.Local().Begin().X()+1, rhs.Local().Begin().Y()+1, rhs.Local().End().Z()-1),
     161                 Index(rhs.Local().End().X()-1, rhs.Local().End().Y()-1, rhs.Local().End().Z()),
     162                 prefactor_inv, diag_inv, off);
    76163}
Note: See TracChangeset for help on using the changeset viewer.