source: src/smoother/gsrb.cpp@ ac6d04

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

Merge recent changes of the vmg library into ScaFaCos.

Includes a fix for the communication problems on Jugene.

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

  • Property mode set to 100644
File size: 4.7 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 <mpi.h>
15
16#include <sstream>
17
18#include "base/discretization.hpp"
19#include "base/stencil.hpp"
20#include "comm/comm.hpp"
21#include "grid/grid.hpp"
22#include "smoother/gsrb.hpp"
23
24using namespace VMG;
25
26static inline void ComputePartial(Grid& sol, Grid& rhs, const Stencil& mat,
27 const Index& begin, const Index& end,
28 const vmg_float& prefactor, const vmg_float& diag_inv,
29 const int& off)
30{
31 int i,j,k;
32 vmg_float temp;
33 Stencil::iterator iter;
34
35 for (i=begin.X(); i<end.X(); ++i)
36 for (j=begin.Y(); j<end.Y(); ++j) {
37
38 int z_begin = begin.Z() + (i + j + begin.Z() + off) % 2;
39
40#ifdef DEBUG
41 int off_sum = MG::GetComm()->LevelSum(rhs, z_begin - begin.Z());
42 assert(z_begin - begin.Z() == 0 || z_begin - begin.Z() == 1);
43 assert(off_sum == 0 || off_sum == MG::GetComm()->Size(rhs));
44#endif /* DEBUG */
45
46 for (k=z_begin; k<end.Z(); k+=2) {
47
48 temp = prefactor * rhs.GetVal(i,j,k);
49
50 for (iter=mat.begin(); iter!=mat.end(); ++iter)
51 temp -= iter->Val() * sol.GetVal(i+iter->Disp().X(),
52 j+iter->Disp().Y(),
53 k+iter->Disp().Z());
54
55 sol(i,j,k) = temp * diag_inv;
56
57 }
58 }
59}
60
61void GaussSeidelRB::Compute(Grid& sol, Grid& rhs)
62{
63 const Stencil& mat = MG::GetDiscretization()->GetStencil();
64 const vmg_float prefactor_inv = 1.0 / MG::GetDiscretization()->OperatorPrefactor(sol);
65 const vmg_float diag_inv = 1.0 / mat.GetDiag();
66 const int off = rhs.Global().LocalBegin().Sum() - rhs.Local().HaloSize1().Sum();
67 const LocalIndices& local = rhs.Local();
68 Comm& comm = *MG::GetComm();
69
70 /*
71 * Compute first halfstep
72 */
73
74 // Start asynchronous communication
75 comm.CommToGhostsAsyncStart(sol);
76
77 // Smooth part not depending on ghost cells
78 ComputePartial(sol, rhs, mat,
79 local.Begin()+1, local.End()-1,
80 prefactor_inv, diag_inv, off+1);
81
82 // Finish asynchronous communication
83 comm.CommToGhostsAsyncFinish(sol);
84
85 /*
86 * Smooth near boundary cells
87 */
88
89 ComputePartial(sol, rhs, mat,
90 local.Begin(),
91 Index(local.Begin().X()+1, local.End().Y(), local.End().Z()),
92 prefactor_inv, diag_inv, off+1);
93
94 ComputePartial(sol, rhs, mat,
95 Index(local.End().X()-1, local.Begin().Y(), local.Begin().Z()),
96 local.End(),
97 prefactor_inv, diag_inv, off+1);
98
99 ComputePartial(sol, rhs, mat,
100 Index(local.Begin().X()+1, local.Begin().Y(), local.Begin().Z()),
101 Index(local.End().X()-1, local.Begin().Y()+1, local.End().Z()),
102 prefactor_inv, diag_inv, off+1);
103
104 ComputePartial(sol, rhs, mat,
105 Index(local.Begin().X()+1, local.End().Y()-1, local.Begin().Z()),
106 Index(local.End().X()-1, local.End().Y(), local.End().Z()),
107 prefactor_inv, diag_inv, off+1);
108
109 ComputePartial(sol, rhs, mat,
110 Index(local.Begin().X()+1, local.Begin().Y()+1, local.Begin().Z()),
111 Index(local.End().X()-1, local.End().Y()-1, local.Begin().Z()+1),
112 prefactor_inv, diag_inv, off+1);
113
114 ComputePartial(sol, rhs, mat,
115 Index(local.Begin().X()+1, local.Begin().Y()+1, local.End().Z()-1),
116 Index(local.End().X()-1, local.End().Y()-1, local.End().Z()),
117 prefactor_inv, diag_inv, off+1);
118
119 /*
120 * Compute second halfstep
121 */
122
123 // Start asynchronous communication
124 comm.CommToGhostsAsyncStart(sol);
125
126 // Smooth part not depending on ghost cells
127 ComputePartial(sol, rhs, mat,
128 local.Begin()+1, local.End()-1,
129 prefactor_inv, diag_inv, off);
130
131 // Finish asynchronous communication
132 comm.CommToGhostsAsyncFinish(sol);
133
134 /*
135 * Smooth near boundary cells
136 */
137
138 ComputePartial(sol, rhs, mat,
139 local.Begin(),
140 Index(local.Begin().X()+1, local.End().Y(), local.End().Z()),
141 prefactor_inv, diag_inv, off);
142
143 ComputePartial(sol, rhs, mat,
144 Index(local.End().X()-1, local.Begin().Y(), local.Begin().Z()),
145 local.End(),
146 prefactor_inv, diag_inv, off);
147
148 ComputePartial(sol, rhs, mat,
149 Index(local.Begin().X()+1, local.Begin().Y(), local.Begin().Z()),
150 Index(local.End().X()-1, local.Begin().Y()+1, local.End().Z()),
151 prefactor_inv, diag_inv, off);
152
153 ComputePartial(sol, rhs, mat,
154 Index(local.Begin().X()+1, local.End().Y()-1, local.Begin().Z()),
155 Index(local.End().X()-1, local.End().Y(), local.End().Z()),
156 prefactor_inv, diag_inv, off);
157
158 ComputePartial(sol, rhs, mat,
159 Index(local.Begin().X()+1, local.Begin().Y()+1, local.Begin().Z()),
160 Index(local.End().X()-1, local.End().Y()-1, local.Begin().Z()+1),
161 prefactor_inv, diag_inv, off);
162
163 ComputePartial(sol, rhs, mat,
164 Index(local.Begin().X()+1, local.Begin().Y()+1, local.End().Z()-1),
165 Index(local.End().X()-1, local.End().Y()-1, local.End().Z()),
166 prefactor_inv, diag_inv, off);
167}
Note: See TracBrowser for help on using the repository browser.