source: src/interface/interface_fcs.cpp@ 51e793

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

vmg: Implement fourth-order discretization of the Poisson equation.

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

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