source: src/interface/interface_fcs.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: 7.4 KB
Line 
1/**
2 * @file interface_fcs.cpp
3 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
4 * @date Mon Apr 18 12:56:20 2011
5 *
6 * @brief VMG::InterfaceFCS
7 *
8 */
9
10#ifdef HAVE_CONFIG_H
11#include <config.h>
12#endif
13
14#ifndef HAVE_MPI
15#error MPI is needed to use the Scafacos interface.
16#endif
17
18#include <mpi.h>
19#ifdef HAVE_MARMOT
20#include <enhancempicalls.h>
21#include <sourceinfompicalls.h>
22#endif
23
24#include "base/object.hpp"
25#include "base/timer.hpp"
26#include "comm/comm_mpi.hpp"
27#include "comm/domain_decomposition_mpi.hpp"
28#include "comm/mpi/error_handler.hpp"
29#include "interface/interface_fcs.h"
30#include "interface/interface_particles.hpp"
31#include "level/level_operator_cs.hpp"
32#include "level/level_operator_fas.hpp"
33#include "samples/discretization_poisson_fd.hpp"
34#include "samples/discretization_poisson_fv.hpp"
35#include "samples/stencils.hpp"
36#include "samples/techniques.hpp"
37#include "smoother/gsrb.cpp"
38#ifdef HAVE_LAPACK
39#include "solver/dgesv.hpp"
40#include "solver/dsysv.hpp"
41#else
42#include "solver/givens.hpp"
43#endif
44#include "solver/solver_regular.hpp"
45#include "solver/solver_singular.hpp"
46
47using namespace VMG;
48
49namespace VMGBackupSettings
50{
51 static vmg_int level = -1;
52 static vmg_int periodic[3] = {-1, -1, -1};
53 static vmg_int max_iter = -1;
54 static vmg_int smoothing_steps = -1;
55 static vmg_int gamma = -1;
56 static vmg_float precision = -1;
57 static vmg_float box_offset[3];
58 static vmg_float box_size = -1.0;
59 static vmg_int near_field_cells = -1;
60 static MPI_Comm mpi_comm;
61}
62
63static void VMG_fcs_init(vmg_int level, vmg_int* periodic,vmg_int max_iter,
64 vmg_int smoothing_steps, vmg_int gamma, vmg_float precision,
65 vmg_float* box_offset, vmg_float box_size,
66 vmg_int near_field_cells, MPI_Comm mpi_comm)
67{
68 VMGBackupSettings::level = level;
69 std::memcpy(VMGBackupSettings::periodic, periodic, 3*sizeof(vmg_int));
70 VMGBackupSettings::max_iter = max_iter;
71 VMGBackupSettings::smoothing_steps = smoothing_steps;
72 VMGBackupSettings::gamma = gamma;
73 VMGBackupSettings::precision = precision;
74 std::memcpy(VMGBackupSettings::box_offset, box_offset, 3*sizeof(vmg_float));
75 VMGBackupSettings::box_size = box_size;
76 VMGBackupSettings::near_field_cells = near_field_cells;
77 VMGBackupSettings::mpi_comm = mpi_comm;
78
79 MPI_Errhandler mpiErrorHandler;
80 MPI_Comm_create_errhandler(VMG::MPI::ConvertToException, &mpiErrorHandler);
81 MPI_Comm_set_errhandler(MPI_COMM_WORLD, mpiErrorHandler);
82
83 const Boundary boundary(periodic[0] ? Periodic : Open,
84 periodic[1] ? Periodic : Open,
85 periodic[2] ? Periodic : Open);
86
87 const bool singular = periodic[0] * periodic[1] * periodic[2];
88
89 /*
90 * Register communication class.
91 */
92 Comm* comm = new CommMPI(boundary, new DomainDecompositionMPI());
93 comm->Register("COMM");
94
95 /*
96 * Register particle interface.
97 */
98 Interface* interface;
99 if (singular)
100 interface = new InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells, 0, 1.0);
101 else
102 interface = new InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells, 2, 1.6);
103 MG::SetInterface(interface, comm);
104
105 /*
106 * Define discretization of the Poisson equation.
107 * Finite volumes for locally refined grids and
108 * finite differences otherwise.
109 */
110 Discretization* discretization;
111 if (singular)
112 discretization = new DiscretizationPoissonFD();
113 else
114 discretization = new DiscretizationPoissonFV();
115 discretization->Register("DISCRETIZATION");
116
117 /*
118 * Use a correction scheme multigrid algorithm with full weight stencils.
119 * Since we use the Gauss-Seidel RB smoother here, this may be replaced
120 * with a half-weight version once debugging is finished.
121 */
122 LevelOperator* lop;
123 if (singular)
124 lop = new LevelOperatorCS(Stencils::RestrictionFullWeight, Stencils::InterpolationTrilinear);
125 else
126 lop = new LevelOperatorFAS(Stencils::RestrictionFullWeight, Stencils::Injection, Stencils::InterpolationTrilinear);
127 lop->Register("LEVEL_OPERATOR");
128
129 /*
130 * Use Gauss-Seidel Red-Black ordering
131 */
132 Smoother* smoother = new GaussSeidelRB();
133 smoother->Register("SMOOTHER");
134
135 /*
136 * Use LAPACK solver when available, otherwise use Givens rotations.
137 */
138#ifdef HAVE_LAPACK
139 Solver* solver = (singular ?
140 static_cast<Solver*>(new DSYSV<SolverSingular>()) :
141 static_cast<Solver*>(new DGESV<SolverRegular>()));
142#else
143 Solver* solver = (singular ?
144 static_cast<Solver*>(new Givens<SolverSingular>()) :
145 static_cast<Solver*>(new Givens<SolverRegular>()));
146#endif
147 solver->Register("SOLVER");
148
149 /*
150 * Set commands for the actual multigrid cycle
151 */
152 if (singular)
153 Techniques::SetCorrectionSchemePeriodic(interface->MinLevel(), interface->MaxLevel(), gamma);
154 else
155 Techniques::SetFullApproximationSchemeDirichlet(interface->MinLevel(), interface->MaxLevel(), gamma);
156
157 /*
158 * Register required parameters
159 */
160 new ObjectStorage<int>("PRESMOOTHSTEPS", smoothing_steps);
161 new ObjectStorage<int>("POSTSMOOTHSTEPS", smoothing_steps);
162 new ObjectStorage<vmg_float>("PRECISION", precision);
163 new ObjectStorage<int>("MAX_ITERATION", max_iter);
164 new ObjectStorage<int>("PARTICLE_NEAR_FIELD_CELLS", near_field_cells);
165
166 /*
167 * Check whether the library is correctly initialized now.
168 */
169 MG::IsInitialized();
170
171 /*
172 * Post init communication class
173 */
174 MG::PostInit();
175}
176
177void VMG_fcs_setup(vmg_int level, vmg_int* periodic, vmg_int max_iter,
178 vmg_int smoothing_steps, vmg_int gamma, vmg_float precision,
179 vmg_float* box_offset, vmg_float box_size,
180 vmg_int near_field_cells, MPI_Comm mpi_comm)
181{
182 if (VMGBackupSettings::level != level ||
183 VMGBackupSettings::periodic[0] != periodic[0] ||
184 VMGBackupSettings::periodic[1] != periodic[1] ||
185 VMGBackupSettings::periodic[2] != periodic[2] ||
186 VMGBackupSettings::max_iter != max_iter ||
187 VMGBackupSettings::smoothing_steps != smoothing_steps ||
188 VMGBackupSettings::gamma != gamma ||
189 VMGBackupSettings::precision != precision ||
190 VMGBackupSettings::box_offset[0] != box_offset[0] ||
191 VMGBackupSettings::box_offset[1] != box_offset[1] ||
192 VMGBackupSettings::box_offset[2] != box_offset[2] ||
193 VMGBackupSettings::box_size != box_size ||
194 VMGBackupSettings::near_field_cells != near_field_cells ||
195 VMGBackupSettings::mpi_comm != mpi_comm) {
196
197 VMG_fcs_destroy();
198 VMG_fcs_init(level, periodic, max_iter,
199 smoothing_steps, gamma, precision,
200 box_offset, box_size, near_field_cells,
201 mpi_comm);
202
203 }
204}
205
206int VMG_fcs_check()
207{
208 const int& near_field_cells = MG::GetFactory().GetObjectStorageVal<int>("PARTICLE_NEAR_FIELD_CELLS");
209 const Multigrid& multigrid = *MG::GetRhs();
210 const Grid& grid = multigrid(multigrid.MaxLevel());
211
212 if (!grid.Global().LocalSize().IsComponentwiseGreater(near_field_cells))
213 return 1;
214
215 return 0;
216}
217
218void VMG_fcs_run(vmg_float* x, vmg_float* q, vmg_float* p, vmg_float* f, vmg_int num_particles_local)
219{
220 /*
221 * Register parameters for later use.
222 */
223 new ObjectStorage<vmg_float*>("PARTICLE_POS_ARRAY", x);
224 new ObjectStorage<vmg_float*>("PARTICLE_CHARGE_ARRAY", q);
225 new ObjectStorage<vmg_float*>("PARTICLE_POTENTIAL_ARRAY", p);
226 new ObjectStorage<vmg_float*>("PARTICLE_FORCE_ARRAY", f);
227 new ObjectStorage<vmg_int>("PARTICLE_NUM_LOCAL", num_particles_local);
228
229 /*
230 * Start the multigrid solver
231 */
232 MG::Solve();
233
234#ifdef DEBUG_MEASURE_TIME
235 Timer::Print();
236#endif
237}
238
239void VMG_fcs_print_timer()
240{
241 Timer::Print();
242}
243
244void VMG_fcs_destroy(void)
245{
246 /*
247 * Delete all data.
248 */
249 MG::Destroy();
250}
Note: See TracBrowser for help on using the repository browser.