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

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

Simplified API.

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

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