source: src/interface/interface_fcs.cpp@ 203547

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

Major vmg update.

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

  • Property mode set to 100644
File size: 6.1 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>
19
20#include "base/object.hpp"
[dfed1c]21#include "base/timer.hpp"
22#include "comm/comm_mpi.hpp"
23#include "comm/domain_decomposition_mpi.hpp"
[48b662]24#include "interface/interface_fcs.h"
25#include "interface/interface_particles.hpp"
26#include "level/level_operator_cs.hpp"
27#include "samples/discretization_poisson_fd.hpp"
28#include "samples/stencils.hpp"
29#include "samples/techniques.hpp"
30#include "smoother/gsrb.cpp"
31#ifdef HAVE_LAPACK
32#include "solver/dgesv.hpp"
33#include "solver/dsysv.hpp"
34#else
35#include "solver/givens.hpp"
36#endif
[dfed1c]37#include "solver/solver_regular.hpp"
38#include "solver/solver_singular.hpp"
[48b662]39
40using namespace VMG;
41
[dfed1c]42namespace VMGBackupSettings
[48b662]43{
[dfed1c]44 static vmg_int level = -1;
45 static vmg_int periodic[3] = {-1, -1, -1};
46 static vmg_int max_iter = -1;
47 static vmg_int smoothing_steps = -1;
48 static vmg_int gamma = -1;
49 static vmg_float precision = -1;
50 static vmg_float box_offset[3];
51 static vmg_float box_size = -1.0;
52 static vmg_int near_field_cells = -1;
53 static MPI_Comm mpi_comm;
54}
[48b662]55
[dfed1c]56static void VMG_fcs_init(vmg_int level, vmg_int* periodic,vmg_int max_iter,
57 vmg_int smoothing_steps, vmg_int gamma, vmg_float precision,
58 vmg_float* box_offset, vmg_float box_size,
59 vmg_int near_field_cells, MPI_Comm mpi_comm)
60{
61 VMGBackupSettings::level = level;
62 std::memcpy(VMGBackupSettings::periodic, periodic, 3*sizeof(vmg_int));
63 VMGBackupSettings::max_iter = max_iter;
64 VMGBackupSettings::smoothing_steps = smoothing_steps;
65 VMGBackupSettings::gamma = gamma;
66 VMGBackupSettings::precision = precision;
67 std::memcpy(VMGBackupSettings::box_offset, box_offset, 3*sizeof(vmg_float));
68 VMGBackupSettings::box_size = box_size;
69 VMGBackupSettings::near_field_cells = near_field_cells;
70 VMGBackupSettings::mpi_comm = mpi_comm;
71
72 const Boundary boundary(periodic[0] ? Periodic : Open,
73 periodic[1] ? Periodic : Open,
74 periodic[2] ? Periodic : Open);
75
76 const bool singular = periodic[0] * periodic[1] * periodic[2];
[48b662]77
78 /*
79 * Register communication class.
80 */
[dfed1c]81 Comm* comm = new CommMPI(boundary, new DomainDecompositionMPI());
[48b662]82 comm->Register("COMM");
83
84 /*
85 * Register particle interface.
86 */
[dfed1c]87 Interface* interface = new VMG::InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells);
[48b662]88 MG::SetInterface(interface, comm);
89
90 /*
91 * Use a finite difference discretization of the Poisson equation
92 */
93 Discretization* discretization = new DiscretizationPoissonFD();
94 discretization->Register("DISCRETIZATION");
95
96 /*
97 * Use a correction scheme multigrid algorithm with full weight stencils.
98 * Since we use the Gauss-Seidel RB smoother here, this may be replaced
99 * with a half-weight version once debugging is finished.
100 */
101 LevelOperator* lop = new LevelOperatorCS(Stencils::RestrictionFullWeight, Stencils::InterpolationTrilinear);
102 lop->Register("LEVEL_OPERATOR");
103
104 /*
105 * Use Gauss-Seidel Red-Black ordering
106 */
107 Smoother* smoother = new GaussSeidelRB();
108 smoother->Register("SMOOTHER");
109
110 /*
111 * Use LAPACK solver when available, otherwise use Givens rotations.
112 */
113#ifdef HAVE_LAPACK
[dfed1c]114 Solver* solver = (singular ?
115 static_cast<Solver*>(new DSYSV<SolverSingular>()) :
116 static_cast<Solver*>(new DGESV<SolverRegular>()));
[48b662]117#else
[dfed1c]118 Solver* solver = (singular ?
119 static_cast<Solver*>(new Givens<SolverSingular>()) :
120 static_cast<Solver*>(new Givens<SolverRegular>()));
[48b662]121#endif
122 solver->Register("SOLVER");
123
124 /*
125 * Set commands for the actual multigrid cycle
126 */
[dfed1c]127 if (singular)
[48b662]128 Techniques::SetCorrectionSchemePeriodic(2, level, gamma);
129 else
130 Techniques::SetCorrectionSchemeDirichlet(2, level, gamma);
131
132 /*
133 * Register required parameters
134 */
135 new ObjectStorage<int>("PRESMOOTHSTEPS", smoothing_steps);
136 new ObjectStorage<int>("POSTSMOOTHSTEPS", smoothing_steps);
137 new ObjectStorage<vmg_float>("PRECISION", precision);
138 new ObjectStorage<int>("MAX_ITERATION", max_iter);
139
[dfed1c]140 new ObjectStorage<int>("PARTICLE_NEAR_FIELD_CELLS", near_field_cells);
141
[48b662]142 /*
143 * Check whether the library is correctly initialized now.
144 */
145 MG::IsInitialized();
146}
147
[dfed1c]148void VMG_fcs_setup(vmg_int level, vmg_int* periodic, vmg_int max_iter,
149 vmg_int smoothing_steps, vmg_int gamma, vmg_float precision,
150 vmg_float* box_offset, vmg_float box_size,
151 vmg_int near_field_cells, MPI_Comm mpi_comm)
152{
153 if (VMGBackupSettings::level != level ||
154 VMGBackupSettings::periodic[0] != periodic[0] ||
155 VMGBackupSettings::periodic[1] != periodic[1] ||
156 VMGBackupSettings::periodic[2] != periodic[2] ||
157 VMGBackupSettings::max_iter != max_iter ||
158 VMGBackupSettings::smoothing_steps != smoothing_steps ||
159 VMGBackupSettings::gamma != gamma ||
160 VMGBackupSettings::precision != precision ||
161 VMGBackupSettings::box_offset[0] != box_offset[0] ||
162 VMGBackupSettings::box_offset[1] != box_offset[1] ||
163 VMGBackupSettings::box_offset[2] != box_offset[2] ||
164 VMGBackupSettings::box_size != box_size ||
165 VMGBackupSettings::near_field_cells != near_field_cells ||
166 VMGBackupSettings::mpi_comm != mpi_comm) {
167
168 VMG_fcs_destroy();
169 VMG_fcs_init(level, periodic, max_iter,
170 smoothing_steps, gamma, precision,
171 box_offset, box_size, near_field_cells,
172 mpi_comm);
173
174 }
175}
176
177void VMG_fcs_run(vmg_float* x, vmg_float* q, vmg_float* p, vmg_float* f, vmg_int num_particles_local)
[48b662]178{
179 /*
180 * Register parameters for later use.
181 */
[dfed1c]182 new ObjectStorage<vmg_float*>("PARTICLE_POS_ARRAY", x);
183 new ObjectStorage<vmg_float*>("PARTICLE_CHARGE_ARRAY", q);
184 new ObjectStorage<vmg_float*>("PARTICLE_POTENTIAL_ARRAY", p);
185 new ObjectStorage<vmg_float*>("PARTICLE_FORCE_ARRAY", f);
186 new ObjectStorage<vmg_int>("PARTICLE_NUM_LOCAL", num_particles_local);
[48b662]187
188 /*
189 * Start the multigrid solver
190 */
191 MG::Solve();
192}
193
[dfed1c]194void VMG_fcs_print_timer()
195{
196 Timer::Print();
197}
198
199void VMG_fcs_destroy(void)
[48b662]200{
201 /*
202 * Delete all data.
203 */
204 MG::Destroy();
205}
Note: See TracBrowser for help on using the repository browser.