source: src/units/particle/interface_fcs.cpp@ 8180d8

Last change on this file since 8180d8 was 8180d8, checked in by Julian Iseringhausen <julian.iseringhausen@…>, 13 years ago

Merge stashed open boundary stuff.

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