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

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

vmg: Added license files and headers (GPLv3).

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

  • Property mode set to 100644
File size: 9.1 KB
Line 
1/*
2 * vmg - a versatile multigrid solver
3 * Copyright (C) 2012 Institute for Numerical Simulation, University of Bonn
4 *
5 * vmg is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * vmg is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19/**
20 * @file interface_fcs.cpp
21 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
22 * @date Mon Apr 18 12:56:20 2011
23 *
24 * @brief VMG::InterfaceFCS
25 *
26 */
27
28#ifdef HAVE_CONFIG_H
29#include <config.h>
30#endif
31
32#ifndef HAVE_MPI
33#error MPI is needed to use the Scafacos interface.
34#endif
35
36#include <mpi.h>
37#ifdef HAVE_MARMOT
38#include <enhancempicalls.h>
39#include <sourceinfompicalls.h>
40#endif
41
42#include "base/object.hpp"
43#include "base/timer.hpp"
44#include "comm/domain_decomposition_mpi.hpp"
45#ifdef DEBUG
46#include "comm/mpi/error_handler.hpp"
47#endif
48#include "level/level_operator_cs.hpp"
49#include "level/level_operator_fas.hpp"
50#include "level/stencils.hpp"
51#include "samples/discretization_poisson_fd.hpp"
52#include "samples/discretization_poisson_fv.hpp"
53#include "samples/techniques.hpp"
54#include "smoother/gsrb.hpp"
55#ifdef HAVE_LAPACK
56#include "solver/dgesv.hpp"
57#include "solver/dsysv.hpp"
58#else
59#include "solver/givens.hpp"
60#endif
61#include "solver/solver_regular.hpp"
62#include "solver/solver_singular.hpp"
63#include "units/particle/comm_mpi_particle.hpp"
64#include "units/particle/interface_fcs.h"
65#include "units/particle/interface_particles.hpp"
66
67
68using namespace VMG;
69
70namespace VMGBackupSettings
71{
72 static vmg_int level = -1;
73 static vmg_int periodic[3] = {-1, -1, -1};
74 static vmg_int max_iter = -1;
75 static vmg_int smoothing_steps = -1;
76 static vmg_int cycle_type = -1;
77 static vmg_float precision = -1;
78 static vmg_float box_offset[3];
79 static vmg_float box_size = -1.0;
80 static vmg_int near_field_cells = -1;
81 static vmg_int interpolation_degree = -1;
82 static vmg_int discretization_order = -1;
83 static MPI_Comm mpi_comm;
84}
85
86static void VMG_fcs_init(vmg_int level, vmg_int* periodic,vmg_int max_iter,
87 vmg_int smoothing_steps, vmg_int cycle_type, vmg_float precision,
88 vmg_float* box_offset, vmg_float box_size,
89 vmg_int near_field_cells, vmg_int interpolation_degree,
90 vmg_int discretization_order, MPI_Comm mpi_comm)
91{
92 VMGBackupSettings::level = level;
93 std::memcpy(VMGBackupSettings::periodic, periodic, 3*sizeof(vmg_int));
94 VMGBackupSettings::max_iter = max_iter;
95 VMGBackupSettings::smoothing_steps = smoothing_steps;
96 VMGBackupSettings::cycle_type = cycle_type;
97 VMGBackupSettings::precision = precision;
98 std::memcpy(VMGBackupSettings::box_offset, box_offset, 3*sizeof(vmg_float));
99 VMGBackupSettings::box_size = box_size;
100 VMGBackupSettings::near_field_cells = near_field_cells;
101 VMGBackupSettings::interpolation_degree = interpolation_degree;
102 VMGBackupSettings::discretization_order = discretization_order;
103 VMGBackupSettings::mpi_comm = mpi_comm;
104
105#ifdef DEBUG
106 MPI_Errhandler mpiErrorHandler;
107 MPI_Comm_create_errhandler(VMG::MPI::ConvertToException, &mpiErrorHandler);
108 MPI_Comm_set_errhandler(MPI_COMM_WORLD, mpiErrorHandler);
109#endif
110
111 const Boundary boundary(periodic[0] ? Periodic : Open,
112 periodic[1] ? Periodic : Open,
113 periodic[2] ? Periodic : Open);
114
115 const bool singular = periodic[0] * periodic[1] * periodic[2];
116
117 /*
118 * Register communication class.
119 */
120 Comm* comm = new Particle::CommMPI(boundary, new DomainDecompositionMPI());
121 comm->Register("COMM");
122
123 /*
124 * Register particle interface.
125 */
126 Interface* interface;
127 if (singular)
128 interface = new InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells, 0, 1.0);
129 else
130 interface = new InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells, 2, 1.6);
131 MG::SetInterface(interface, comm);
132
133 /*
134 * Define discretization of the Poisson equation.
135 * Finite volumes for locally refined grids and
136 * finite differences otherwise.
137 */
138 Discretization* discretization;
139 if (singular)
140 discretization = new DiscretizationPoissonFD(discretization_order);
141 else
142 discretization = new DiscretizationPoissonFV(discretization_order);
143 discretization->Register("DISCRETIZATION");
144
145 /*
146 * Use a correction scheme multigrid algorithm with full weight stencils.
147 * Since we use the Gauss-Seidel RB smoother here, this may be replaced
148 * with a half-weight version once debugging is finished.
149 */
150 LevelOperator* lop;
151 if (singular)
152 lop = new LevelOperatorCS(Stencils::RestrictionFullWeight, Stencils::InterpolationTrilinear);
153 else
154 lop = new LevelOperatorFAS(Stencils::RestrictionFullWeight, Stencils::Injection, Stencils::InterpolationTrilinear);
155 lop->Register("LEVEL_OPERATOR");
156
157 /*
158 * Use Gauss-Seidel Red-Black ordering
159 */
160 Smoother* smoother = new GaussSeidelRB();
161 smoother->Register("SMOOTHER");
162
163 /*
164 * Use LAPACK solver when available, otherwise use Givens rotations.
165 */
166#ifdef HAVE_LAPACK
167 Solver* solver = (singular ?
168 static_cast<Solver*>(new DSYSV<SolverSingular>()) :
169 static_cast<Solver*>(new DGESV<SolverRegular>()));
170#else
171 Solver* solver = (singular ?
172 static_cast<Solver*>(new Givens<SolverSingular>()) :
173 static_cast<Solver*>(new Givens<SolverRegular>()));
174#endif
175 solver->Register("SOLVER");
176
177 /*
178 * Set commands for the actual multigrid cycle
179 */
180 if (singular)
181 Techniques::SetCorrectionSchemePeriodic(interface->MinLevel(), interface->MaxLevel(), cycle_type);
182 else
183 Techniques::SetFullApproximationSchemeDirichlet(interface->MinLevel(), interface->MaxLevel(), cycle_type);
184
185 /*
186 * Register required parameters
187 */
188 new ObjectStorage<int>("PRESMOOTHSTEPS", smoothing_steps);
189 new ObjectStorage<int>("POSTSMOOTHSTEPS", smoothing_steps);
190 new ObjectStorage<vmg_float>("PRECISION", precision);
191 new ObjectStorage<int>("MAX_ITERATION", max_iter);
192 new ObjectStorage<int>("PARTICLE_NEAR_FIELD_CELLS", near_field_cells);
193 new ObjectStorage<int>("PARTICLE_INTERPOLATION_DEGREE", interpolation_degree);
194
195 /*
196 * Check whether the library is correctly initialized now.
197 */
198 MG::IsInitialized();
199
200 /*
201 * Post init communication class
202 */
203 MG::PostInit();
204}
205
206void VMG_fcs_setup(vmg_int level, vmg_int* periodic, vmg_int max_iter,
207 vmg_int smoothing_steps, vmg_int cycle_type, vmg_float precision,
208 vmg_float* box_offset, vmg_float box_size,
209 vmg_int near_field_cells, vmg_int interpolation_degree,
210 vmg_int discretization_order, MPI_Comm mpi_comm)
211{
212 if (VMGBackupSettings::level != level ||
213 VMGBackupSettings::periodic[0] != periodic[0] ||
214 VMGBackupSettings::periodic[1] != periodic[1] ||
215 VMGBackupSettings::periodic[2] != periodic[2] ||
216 VMGBackupSettings::max_iter != max_iter ||
217 VMGBackupSettings::smoothing_steps != smoothing_steps ||
218 VMGBackupSettings::cycle_type != cycle_type ||
219 VMGBackupSettings::precision != precision ||
220 VMGBackupSettings::box_offset[0] != box_offset[0] ||
221 VMGBackupSettings::box_offset[1] != box_offset[1] ||
222 VMGBackupSettings::box_offset[2] != box_offset[2] ||
223 VMGBackupSettings::box_size != box_size ||
224 VMGBackupSettings::near_field_cells != near_field_cells ||
225 VMGBackupSettings::interpolation_degree != interpolation_degree ||
226 VMGBackupSettings::discretization_order != discretization_order ||
227 VMGBackupSettings::mpi_comm != mpi_comm) {
228
229 VMG_fcs_destroy();
230 VMG_fcs_init(level, periodic, max_iter,
231 smoothing_steps, cycle_type, precision,
232 box_offset, box_size, near_field_cells,
233 interpolation_degree, discretization_order,
234 mpi_comm);
235
236 }
237}
238
239int VMG_fcs_check()
240{
241 const int& near_field_cells = MG::GetFactory().GetObjectStorageVal<int>("PARTICLE_NEAR_FIELD_CELLS");
242 const Multigrid& multigrid = *MG::GetRhs();
243 const Grid& grid = multigrid(multigrid.MaxLevel());
244
245 int error_code;
246
247 if (!grid.Global().LocalSize().IsComponentwiseGreater(near_field_cells))
248 error_code = 1;
249 else
250 error_code = 0;
251
252 error_code = MG::GetComm()->GlobalMax(error_code);
253
254 return error_code;
255}
256
257void VMG_fcs_run(vmg_float* x, vmg_float* q, vmg_float* p, vmg_float* f, vmg_int num_particles_local)
258{
259 /*
260 * Register parameters for later use.
261 */
262 new ObjectStorage<vmg_float*>("PARTICLE_POS_ARRAY", x);
263 new ObjectStorage<vmg_float*>("PARTICLE_CHARGE_ARRAY", q);
264 new ObjectStorage<vmg_float*>("PARTICLE_POTENTIAL_ARRAY", p);
265 new ObjectStorage<vmg_float*>("PARTICLE_FIELD_ARRAY", f);
266 new ObjectStorage<vmg_int>("PARTICLE_NUM_LOCAL", num_particles_local);
267
268 /*
269 * Start the multigrid solver
270 */
271 MG::Solve();
272
273#ifdef DEBUG_MEASURE_TIME
274 Timer::Print();
275#endif
276}
277
278void VMG_fcs_print_timer()
279{
280 Timer::Print();
281}
282
283void VMG_fcs_destroy(void)
284{
285 /*
286 * Delete all data.
287 */
288 MG::Destroy();
289}
Note: See TracBrowser for help on using the repository browser.