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

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

vmg: Work on output verbosity.

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

  • Property mode set to 100644
File size: 7.8 KB
RevLine 
[fcf7f6]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
[48b662]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>
[ac6d04]37#ifdef HAVE_MARMOT
38#include <enhancempicalls.h>
39#include <sourceinfompicalls.h>
40#endif
[48b662]41
42#include "base/object.hpp"
[dfed1c]43#include "base/timer.hpp"
44#include "comm/domain_decomposition_mpi.hpp"
[f003a9]45#ifdef DEBUG
[ac6d04]46#include "comm/mpi/error_handler.hpp"
[f003a9]47#endif
[b57b9b]48#include "cycles/cycle_cs_periodic.hpp"
49#include "cycles/cycle_fas_dirichlet.hpp"
[48b662]50#include "level/level_operator_cs.hpp"
[ac6d04]51#include "level/level_operator_fas.hpp"
[b7b317]52#include "level/stencils.hpp"
[48b662]53#include "samples/discretization_poisson_fd.hpp"
[ac6d04]54#include "samples/discretization_poisson_fv.hpp"
[952298]55#include "smoother/gsrb_poisson_2.hpp"
56#include "smoother/gsrb_poisson_4.hpp"
[48b662]57#include "solver/givens.hpp"
[dfed1c]58#include "solver/solver_regular.hpp"
59#include "solver/solver_singular.hpp"
[f003a9]60#include "units/particle/comm_mpi_particle.hpp"
61#include "units/particle/interface_fcs.h"
62#include "units/particle/interface_particles.hpp"
63
[48b662]64
65using namespace VMG;
66
[dfed1c]67namespace VMGBackupSettings
[48b662]68{
[dfed1c]69 static vmg_int level = -1;
70 static vmg_int periodic[3] = {-1, -1, -1};
71 static vmg_int max_iter = -1;
72 static vmg_int smoothing_steps = -1;
[728149]73 static vmg_int cycle_type = -1;
[dfed1c]74 static vmg_float precision = -1;
75 static vmg_float box_offset[3];
76 static vmg_float box_size = -1.0;
77 static vmg_int near_field_cells = -1;
[a150d0]78 static vmg_int interpolation_degree = -1;
[4571da]79 static vmg_int discretization_order = -1;
[dfed1c]80 static MPI_Comm mpi_comm;
81}
[48b662]82
[dfed1c]83static void VMG_fcs_init(vmg_int level, vmg_int* periodic,vmg_int max_iter,
[728149]84 vmg_int smoothing_steps, vmg_int cycle_type, vmg_float precision,
[dfed1c]85 vmg_float* box_offset, vmg_float box_size,
[a150d0]86 vmg_int near_field_cells, vmg_int interpolation_degree,
[4571da]87 vmg_int discretization_order, MPI_Comm mpi_comm)
[dfed1c]88{
89 VMGBackupSettings::level = level;
90 std::memcpy(VMGBackupSettings::periodic, periodic, 3*sizeof(vmg_int));
91 VMGBackupSettings::max_iter = max_iter;
92 VMGBackupSettings::smoothing_steps = smoothing_steps;
[728149]93 VMGBackupSettings::cycle_type = cycle_type;
[dfed1c]94 VMGBackupSettings::precision = precision;
95 std::memcpy(VMGBackupSettings::box_offset, box_offset, 3*sizeof(vmg_float));
96 VMGBackupSettings::box_size = box_size;
97 VMGBackupSettings::near_field_cells = near_field_cells;
[a150d0]98 VMGBackupSettings::interpolation_degree = interpolation_degree;
[4571da]99 VMGBackupSettings::discretization_order = discretization_order;
[dfed1c]100 VMGBackupSettings::mpi_comm = mpi_comm;
101
[379700]102#ifdef DEBUG
[ac6d04]103 MPI_Errhandler mpiErrorHandler;
104 MPI_Comm_create_errhandler(VMG::MPI::ConvertToException, &mpiErrorHandler);
[2a5451]105 MPI_Comm_set_errhandler(mpi_comm, mpiErrorHandler);
[379700]106#endif
[ac6d04]107
[dfed1c]108 const Boundary boundary(periodic[0] ? Periodic : Open,
109 periodic[1] ? Periodic : Open,
110 periodic[2] ? Periodic : Open);
111
112 const bool singular = periodic[0] * periodic[1] * periodic[2];
[48b662]113
114 /*
[a24b80]115 * Choose multigrid components
[48b662]116 */
[a24b80]117 if (singular) {
[48b662]118
[b57b9b]119 new Particle::CommMPI(boundary, new DomainDecompositionMPI(), mpi_comm);
120 new DiscretizationPoissonFD(discretization_order);
121 new InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells, 0, 1.0);
122 new LevelOperatorCS(Stencils::RestrictionFullWeight, Stencils::InterpolationTrilinear);
123 new Givens<SolverSingular>();
124 new CycleCSPeriodic(cycle_type);
[48b662]125
[952298]126 }else {
[a24b80]127
[b57b9b]128 new Particle::CommMPI(boundary, new DomainDecompositionMPI(), mpi_comm);
129 new DiscretizationPoissonFV(discretization_order);
130 new InterfaceParticles(boundary, 2, level, Vector(box_offset), box_size, near_field_cells, 2, 1.6);
131 new LevelOperatorFAS(Stencils::RestrictionFullWeight, Stencils::Injection, Stencils::InterpolationTrilinear);
132 new Givens<SolverRegular>();
133 new CycleFASDirichlet(cycle_type);
[a24b80]134
[952298]135 }
[48b662]136
137 /*
[a24b80]138 * Use Gauss-Seidel Red-Black ordering
[48b662]139 */
[a24b80]140 if (discretization_order == 2)
[b57b9b]141 new GaussSeidelRBPoisson2();
[a24b80]142 else
[b57b9b]143 new GaussSeidelRBPoisson4();
[48b662]144
145 /*
146 * Register required parameters
147 */
148 new ObjectStorage<int>("PRESMOOTHSTEPS", smoothing_steps);
149 new ObjectStorage<int>("POSTSMOOTHSTEPS", smoothing_steps);
150 new ObjectStorage<vmg_float>("PRECISION", precision);
151 new ObjectStorage<int>("MAX_ITERATION", max_iter);
[dfed1c]152 new ObjectStorage<int>("PARTICLE_NEAR_FIELD_CELLS", near_field_cells);
[a150d0]153 new ObjectStorage<int>("PARTICLE_INTERPOLATION_DEGREE", interpolation_degree);
[dfed1c]154
[48b662]155 /*
[a24b80]156 * Post init
[48b662]157 */
[a24b80]158 MG::PostInit();
[894a5f]159
160 /*
[a24b80]161 * Check whether the library is correctly initialized now.
[894a5f]162 */
[a24b80]163 MG::IsInitialized();
[48b662]164}
165
[dfed1c]166void VMG_fcs_setup(vmg_int level, vmg_int* periodic, vmg_int max_iter,
[728149]167 vmg_int smoothing_steps, vmg_int cycle_type, vmg_float precision,
[dfed1c]168 vmg_float* box_offset, vmg_float box_size,
[a150d0]169 vmg_int near_field_cells, vmg_int interpolation_degree,
[4571da]170 vmg_int discretization_order, MPI_Comm mpi_comm)
[dfed1c]171{
172 if (VMGBackupSettings::level != level ||
173 VMGBackupSettings::periodic[0] != periodic[0] ||
174 VMGBackupSettings::periodic[1] != periodic[1] ||
175 VMGBackupSettings::periodic[2] != periodic[2] ||
176 VMGBackupSettings::max_iter != max_iter ||
177 VMGBackupSettings::smoothing_steps != smoothing_steps ||
[728149]178 VMGBackupSettings::cycle_type != cycle_type ||
[dfed1c]179 VMGBackupSettings::precision != precision ||
180 VMGBackupSettings::box_offset[0] != box_offset[0] ||
181 VMGBackupSettings::box_offset[1] != box_offset[1] ||
182 VMGBackupSettings::box_offset[2] != box_offset[2] ||
183 VMGBackupSettings::box_size != box_size ||
184 VMGBackupSettings::near_field_cells != near_field_cells ||
[a150d0]185 VMGBackupSettings::interpolation_degree != interpolation_degree ||
[4571da]186 VMGBackupSettings::discretization_order != discretization_order ||
[dfed1c]187 VMGBackupSettings::mpi_comm != mpi_comm) {
188
189 VMG_fcs_destroy();
190 VMG_fcs_init(level, periodic, max_iter,
[728149]191 smoothing_steps, cycle_type, precision,
[dfed1c]192 box_offset, box_size, near_field_cells,
[4571da]193 interpolation_degree, discretization_order,
194 mpi_comm);
[dfed1c]195
196 }
197}
198
[ac6d04]199int VMG_fcs_check()
200{
201 const int& near_field_cells = MG::GetFactory().GetObjectStorageVal<int>("PARTICLE_NEAR_FIELD_CELLS");
202 const Multigrid& multigrid = *MG::GetRhs();
203 const Grid& grid = multigrid(multigrid.MaxLevel());
204
[e3da7f]205 vmg_int error_code = 0;
[7dd744]206
[ac6d04]207 if (!grid.Global().LocalSize().IsComponentwiseGreater(near_field_cells))
[7dd744]208 error_code = 1;
[ac6d04]209
[ef94e7]210 return MG::GetComm()->GlobalMax(error_code);
[ac6d04]211}
212
[dfed1c]213void VMG_fcs_run(vmg_float* x, vmg_float* q, vmg_float* p, vmg_float* f, vmg_int num_particles_local)
[48b662]214{
215 /*
216 * Register parameters for later use.
217 */
[dfed1c]218 new ObjectStorage<vmg_float*>("PARTICLE_POS_ARRAY", x);
219 new ObjectStorage<vmg_float*>("PARTICLE_CHARGE_ARRAY", q);
220 new ObjectStorage<vmg_float*>("PARTICLE_POTENTIAL_ARRAY", p);
[716da7]221 new ObjectStorage<vmg_float*>("PARTICLE_FIELD_ARRAY", f);
[dfed1c]222 new ObjectStorage<vmg_int>("PARTICLE_NUM_LOCAL", num_particles_local);
[48b662]223
224 /*
225 * Start the multigrid solver
226 */
227 MG::Solve();
[894a5f]228
229 Timer::Print();
[48b662]230}
231
[dfed1c]232void VMG_fcs_print_timer()
233{
234 Timer::Print();
235}
236
237void VMG_fcs_destroy(void)
[48b662]238{
239 /*
240 * Delete all data.
241 */
242 MG::Destroy();
243}
Note: See TracBrowser for help on using the repository browser.