source: src/units/particle/interface_fcs.cpp@ b7b317

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

Fixed newly introduced error due to unexpected behaviour of git-svn. Sorry.

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

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