source: src/interface/interface_fcs.cpp@ 728149

Last change on this file since 728149 was 728149, checked in by Frederik Heber <heber@…>, 14 years ago

Conformalized vmg getter functions with other solvers.

  • also rename gamma -> cycle_type.
  • all getters now return FCSResult and require two parameters.

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

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