source: src/interface/interface_fcs.cpp@ 71b148

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

vmg: When internally checking for configuration errors, return the maximum error code of all processes on all processes instead of letting each process return it's own error code.

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

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