source: src/interface/interface_fcs.cpp@ 48b662

Last change on this file since 48b662 was 48b662, checked in by Olaf Lenz <olenz@…>, 14 years ago

Moved files in scafacos_fcs one level up.

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

  • Property mode set to 100644
File size: 4.3 KB
Line 
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"
21#include "comm/comm_serial.hpp"
22#include "interface/interface_fcs.h"
23#include "interface/interface_particles.hpp"
24#include "level/level_operator_cs.hpp"
25#include "samples/discretization_poisson_fd.hpp"
26#include "samples/stencils.hpp"
27#include "samples/techniques.hpp"
28#include "smoother/gsrb.cpp"
29#ifdef HAVE_LAPACK
30#include "solver/dgesv.hpp"
31#include "solver/dsysv.hpp"
32#else
33#include "solver/givens.hpp"
34#endif
35#include "solver/solver_dirichlet.hpp"
36#include "solver/solver_periodic.hpp"
37
38using namespace VMG;
39
40void VMGInterfaceFCSInit(vmg_int level, vmg_int periodic, vmg_int spline_degree,
41 vmg_int max_iter, vmg_int smoothing_steps, vmg_int gamma,
42 vmg_float precision, vmg_float box_begin, vmg_float box_end,
43 MPI_Comm mpi_comm)
44{
45 const BC bc = (periodic ? Periodic : Dirichlet);
46
47 /*
48 * Register all the parameters for later use.
49 */
50 new ObjectStorage<int>("LEVEL_FCS", level);
51 new ObjectStorage<int>("PERIODIC_FCS", periodic);
52 new ObjectStorage<int>("SPLINE_DEGREE_FCS", spline_degree);
53 new ObjectStorage<int>("MAX_ITER_FCS", max_iter);
54 new ObjectStorage<int>("SMOOTHING_STEPS_FCS", smoothing_steps);
55 new ObjectStorage<int>("GAMMA_FCS", gamma);
56 new ObjectStorage<vmg_float>("PRECISION_FCS", precision);
57 new ObjectStorage<vmg_float>("BOX_BEGIN_FCS", box_begin);
58 new ObjectStorage<vmg_float>("BOX_END_FCS", box_end);
59 new ObjectStorage<MPI_Comm>("MPI_COMM_FCS", mpi_comm);
60
61 /*
62 * Register communication class.
63 * For now, parallel version is not yet released.
64 */
65 Comm* comm = new CommSerial(bc);
66 comm->Register("COMM");
67
68 /*
69 * Register particle interface.
70 */
71 Interface* interface = new VMG::InterfaceParticles(bc, 2, level, box_begin, box_end);
72 MG::SetInterface(interface, comm);
73
74 /*
75 * Use a finite difference discretization of the Poisson equation
76 */
77 Discretization* discretization = new DiscretizationPoissonFD();
78 discretization->Register("DISCRETIZATION");
79
80 /*
81 * Use a correction scheme multigrid algorithm with full weight stencils.
82 * Since we use the Gauss-Seidel RB smoother here, this may be replaced
83 * with a half-weight version once debugging is finished.
84 */
85 LevelOperator* lop = new LevelOperatorCS(Stencils::RestrictionFullWeight, Stencils::InterpolationTrilinear);
86 lop->Register("LEVEL_OPERATOR");
87
88 /*
89 * Use Gauss-Seidel Red-Black ordering
90 */
91 Smoother* smoother = new GaussSeidelRB();
92 smoother->Register("SMOOTHER");
93
94 /*
95 * Use LAPACK solver when available, otherwise use Givens rotations.
96 */
97#ifdef HAVE_LAPACK
98 Solver* solver = (periodic ?
99 static_cast<Solver*>(new DSYSV<SolverPeriodic>()) :
100 static_cast<Solver*>(new DGESV<SolverDirichlet>()));
101#else
102 Solver* solver = (periodic ?
103 static_cast<Solver*>(new Givens<SolverPeriodic>()) :
104 static_cast<Solver*>(new Givens<SolverDirichlet>()));
105#endif
106 solver->Register("SOLVER");
107
108 /*
109 * Set commands for the actual multigrid cycle
110 */
111 if (periodic)
112 Techniques::SetCorrectionSchemePeriodic(2, level, gamma);
113 else
114 Techniques::SetCorrectionSchemeDirichlet(2, level, gamma);
115
116 /*
117 * Register required parameters
118 */
119 new ObjectStorage<int>("PRESMOOTHSTEPS", smoothing_steps);
120 new ObjectStorage<int>("POSTSMOOTHSTEPS", smoothing_steps);
121 new ObjectStorage<vmg_float>("PRECISION", precision);
122 new ObjectStorage<int>("MAX_ITERATION", max_iter);
123
124 /*
125 * Check whether the library is correctly initialized now.
126 */
127 MG::IsInitialized();
128}
129
130void VMGInterfaceFCSRun(vmg_float* x, vmg_float* q, vmg_float* p, vmg_float* f, vmg_int num_particles_local)
131{
132 /*
133 * Register parameters for later use.
134 */
135 new ObjectStorage<vmg_float*>("X_FCS", x);
136 new ObjectStorage<vmg_float*>("Q_FCS", q);
137 new ObjectStorage<vmg_float*>("P_FCS", p);
138 new ObjectStorage<vmg_float*>("F_FCS", f);
139 new ObjectStorage<vmg_int>("NUM_PARTICLES_LOCAL_FCS", num_particles_local);
140
141 /*
142 * Start the multigrid solver
143 */
144 MG::Solve();
145}
146
147void VMGInterfaceFCSDestroy(void)
148{
149 /*
150 * Delete all data.
151 */
152 MG::Destroy();
153}
Note: See TracBrowser for help on using the repository browser.