source: src/smoother/gsrb.cpp@ 894a5f

Last change on this file since 894a5f was 894a5f, checked in by Julian Iseringhausen <isering@…>, 14 years ago

Parallel performance update.

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

  • Property mode set to 100644
File size: 4.8 KB
Line 
1/**
2 * @file gsrb.cpp
3 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
4 * @date Mon Apr 18 13:08:20 2011
5 *
6 * @brief Gauss-Seidel Red Black method
7 *
8 */
9
10#ifdef HAVE_CONFIG_H
11#include <config.h>
12#endif
13
14#include <sstream>
15
16#include "base/discretization.hpp"
17#include "base/stencil.hpp"
18#include "comm/comm.hpp"
19#include "grid/grid.hpp"
20#include "smoother/gsrb.hpp"
21
22using 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}
48
49void GaussSeidelRB::Compute(Grid& sol, Grid& rhs)
50{
51#ifdef DEBUG_MATRIX_CHECKS
52 sol.IsConsistent();
53 rhs.IsConsistent();
54#endif
55
56#ifdef DEBUG
57 sol.IsCompatible(rhs);
58#endif
59
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();
65
66 /*
67 * Compute first halfstep
68 */
69
70 // Start asynchronous communication
71 comm.CommToGhostsAsyncStart(sol);
72
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);
77
78 // Finish asynchronous communication
79 comm.CommToGhostsAsyncFinish();
80
81 /*
82 * Smooth near boundary cells
83 */
84
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);
89
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);
94
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);
99
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);
104
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);
109
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);
114
115 /*
116 * Compute second halfstep
117 */
118
119 // Start asynchronous communication
120 comm.CommToGhostsAsyncStart(sol);
121
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);
126
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);
163}
Note: See TracBrowser for help on using the repository browser.