source: src/mg.cpp@ f003a9

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

Refactored vmg in order to separate the core library and the particle simulation part properly.

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

  • Property mode set to 100644
File size: 6.7 KB
Line 
1/**
2 * @file mg.cpp
3 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
4 * @date Sat Jun 12 20:36:24 2010
5 *
6 * @brief A multigrid solver
7 *
8 * This file contains the implementation of the main class for
9 * a multigrid solver.
10 *
11 */
12
13#ifdef HAVE_CONFIG_H
14#include <config.h>
15#endif
16
17#ifdef DEBUG_MEASURE_TIME
18#ifdef HAVE_MPI
19#include <mpi.h>
20#ifdef HAVE_MARMOT
21#include <enhancempicalls.h>
22#include <sourceinfompicalls.h>
23#endif
24#endif
25#endif
26
27#include <cmath>
28#include <cstdio>
29#include <cfloat>
30#include <ctime>
31#include <string>
32#include <fstream>
33#include <sstream>
34
35#include "base/command_list.hpp"
36#include "base/discretization.hpp"
37#include "base/factory.hpp"
38#include "base/interface.hpp"
39#include "base/timer.hpp"
40#include "comm/comm.hpp"
41#include "grid/grid.hpp"
42#include "grid/tempgrid.hpp"
43#include "level/level_operator.hpp"
44#include "smoother/smoother.hpp"
45#include "solver/solver.hpp"
46#include "mg.hpp"
47
48using namespace VMG;
49
50#define REGISTER_COMMAND(a) extern void Initialize##a();Initialize##a();
51
52VMG::CommandFactory MG::command_factory;
53
54static void VMGRegisterBuiltinCommands()
55{
56 REGISTER_COMMAND(VMGCommandCheckConsistency);
57 REGISTER_COMMAND(VMGCommandCheckIterationCounter);
58 REGISTER_COMMAND(VMGCommandCheckRelativeResidual);
59 REGISTER_COMMAND(VMGCommandCheckResidual);
60 REGISTER_COMMAND(VMGCommandClearCoarseLevels);
61 REGISTER_COMMAND(VMGCommandClearGrid);
62 REGISTER_COMMAND(VMGCommandComputeResidualNorm);
63 REGISTER_COMMAND(VMGCommandCopyBoundary);
64 REGISTER_COMMAND(VMGCommandExecuteCycle);
65 REGISTER_COMMAND(VMGCommandExecuteCycleLoop);
66 REGISTER_COMMAND(VMGCommandExecuteFullCycle);
67 REGISTER_COMMAND(VMGCommandExecuteFullCycleLoop);
68 REGISTER_COMMAND(VMGCommandExportSolution);
69 REGISTER_COMMAND(VMGCommandForceDiscreteCompatibility);
70 REGISTER_COMMAND(VMGCommandImportRightHandSide);
71 REGISTER_COMMAND(VMGCommandInterpolateFMG);
72 REGISTER_COMMAND(VMGCommandInitializeIterationCounter);
73 REGISTER_COMMAND(VMGCommandInitializeResidualNorm);
74 REGISTER_COMMAND(VMGCommandNOP);
75 REGISTER_COMMAND(VMGCommandPrintAllSettings);
76 REGISTER_COMMAND(VMGCommandPrintDefect);
77 REGISTER_COMMAND(VMGCommandPrintGridStructure);
78 REGISTER_COMMAND(VMGCommandPrintGrid);
79 REGISTER_COMMAND(VMGCommandPrintResidualNorm);
80 REGISTER_COMMAND(VMGCommandPrintRunningTime);
81 REGISTER_COMMAND(VMGCommandProlongate);
82 REGISTER_COMMAND(VMGCommandRestrict);
83 REGISTER_COMMAND(VMGCommandSetAverageToZero);
84 REGISTER_COMMAND(VMGCommandSetCoarserDirichletValues);
85 REGISTER_COMMAND(VMGCommandSetLevel);
86 REGISTER_COMMAND(VMGCommandSmooth);
87 REGISTER_COMMAND(VMGCommandSolve);
88}
89
90MG::MG()
91{
92 state = 0;
93 VMGRegisterBuiltinCommands();
94}
95
96MG::~MG()
97{
98 MG::Destroy();
99}
100
101// Brings Multigrid back to starting state.
102void MG::Destroy()
103{
104 MG::Instance()->factories.clear();
105 MG::Instance()->state = 0;
106 Timer::Clear();
107}
108
109void MG::SetInterface(VMG::Interface* interface, VMG::Comm* comm)
110{
111 interface->Register("INTERFACE");
112
113 Multigrid* sol = new Multigrid(comm, interface);
114 Multigrid* rhs = new Multigrid(comm, interface);
115 TempGrid* temp = new TempGrid();
116
117 sol->Register("SOL");
118 rhs->Register("RHS");
119 temp->Register("TEMPGRID");
120
121 new ObjectStorage<int>("GLOBAL_MAXLEVEL", sol->GlobalMaxLevel());
122 new ObjectStorage<int>("MINLEVEL", sol->MinLevel());
123 new ObjectStorage<int>("MAXLEVEL", sol->MaxLevel());
124}
125
126/*
127 * Post init communication class
128 */
129void MG::PostInit()
130{
131 GetComm()->PostInit(*GetSol(), *GetRhs());
132}
133
134/**
135 * Solves a given system with a multigrid method
136 *
137 */
138void MG::Solve()
139{
140#ifdef DEBUG_MEASURE_TIME
141#ifdef HAVE_MPI
142 MPI_Barrier(MPI_COMM_WORLD);
143#endif
144 Timer::Start("CompleteRunningTime");
145#endif
146
147 CommandList* cl_init = MG::GetFactory().Get("COMMANDLIST_INIT")->Cast<CommandList>();
148 CommandList* cl_loop = MG::GetFactory().Get("COMMANDLIST_LOOP")->Cast<CommandList>();
149 CommandList* cl_finalize = MG::GetFactory().Get("COMMANDLIST_FINALIZE")->Cast<CommandList>();
150
151 cl_init->ExecuteList();
152
153 while (cl_loop->ExecuteList() == Continue);
154
155 cl_finalize->ExecuteList();
156
157#ifdef DEBUG_MEASURE_TIME
158#ifdef HAVE_MPI
159 MPI_Barrier(MPI_COMM_WORLD);
160#endif
161 Timer::Stop("CompleteRunningTime");
162#ifdef DEBUG_MEASURE_TIME_OUTPUT
163#ifdef HAVE_MPI
164 Timer::PrintGlobal();
165#else
166 Timer::Print();
167#endif
168#endif
169#endif
170}
171
172void MG::SetState(const int& state_)
173{
174 MG::Instance()->state = state_;
175}
176
177VMG::Factory& MG::GetFactory()
178{
179 std::map<int, VMG::Factory>::iterator iter = MG::Instance()->factories.find(MG::Instance()->state);
180
181 if (iter == MG::Instance()->factories.end())
182 iter = MG::Instance()->factories.insert(std::make_pair(MG::Instance()->state, Factory())).first;
183
184 assert(iter != MG::Instance()->factories.end());
185
186 return iter->second;
187}
188
189VMG::CommandFactory& MG::GetCommands()
190{
191 return MG::command_factory;
192}
193
194Comm* MG::GetComm()
195{
196 return MG::GetFactory().Get("COMM")->Cast<VMG::Comm>();
197}
198
199Discretization* MG::GetDiscretization()
200{
201 return MG::GetFactory().Get("DISCRETIZATION")->Cast<VMG::Discretization>();
202}
203
204LevelOperator* MG::GetLevelOperator()
205{
206 return MG::GetFactory().Get("LEVEL_OPERATOR")->Cast<VMG::LevelOperator>();
207}
208
209Multigrid* MG::GetRhs()
210{
211 return MG::GetFactory().Get("RHS")->Cast<VMG::Multigrid>();
212}
213
214Multigrid* MG::GetSol()
215{
216 return MG::GetFactory().Get("SOL")->Cast<VMG::Multigrid>();
217}
218
219VMG::Grid& MG::GetRhsMaxLevel()
220{
221 return (*MG::GetRhs())(MG::GetRhs()->MaxLevel());
222}
223
224VMG::Grid& MG::GetSolMaxLevel()
225{
226 return (*MG::GetSol())(MG::GetSol()->MaxLevel());
227}
228
229Smoother* MG::GetSmoother()
230{
231 return MG::GetFactory().Get("SMOOTHER")->Cast<VMG::Smoother>();
232}
233
234Solver* MG::GetSolver()
235{
236 return MG::GetFactory().Get("SOLVER")->Cast<VMG::Solver>();
237}
238
239TempGrid* MG::GetTempGrid()
240{
241 return MG::GetFactory().Get("TEMPGRID")->Cast<VMG::TempGrid>();
242}
243
244Interface* MG::GetInterface()
245{
246 return MG::GetFactory().Get("INTERFACE")->Cast<VMG::Interface>();
247}
248
249static bool CheckObject(std::string id)
250{
251 Object *obj = MG::GetFactory().Get(id);
252
253#ifdef DEBUG_OUTPUT
254 if (obj == NULL)
255 printf("\nMultigrid: CLASS %s NOT INITIALIZED\n\n", id.c_str());
256#endif
257
258 return obj != NULL;
259}
260
261bool MG::IsInitialized()
262{
263 bool init = true;
264
265 init &= CheckObject("COMM");
266 init &= CheckObject("LEVEL_OPERATOR");
267 init &= CheckObject("RHS");
268 init &= CheckObject("SOL");
269 init &= CheckObject("SOLVER");
270 init &= CheckObject("SMOOTHER");
271 init &= CheckObject("DISCRETIZATION");
272 init &= CheckObject("MAX_ITERATION");
273 init &= CheckObject("PRECISION");
274 init &= CheckObject("PRESMOOTHSTEPS");
275 init &= CheckObject("POSTSMOOTHSTEPS");
276 init &= CheckObject("COMMANDLIST_INIT");
277 init &= CheckObject("COMMANDLIST_LOOP");
278 init &= CheckObject("COMMANDLIST_FINALIZE");
279 init &= CheckObject("MINLEVEL");
280 init &= CheckObject("MAXLEVEL");
281 init &= CheckObject("GLOBAL_MAXLEVEL");
282 init &= CheckObject("INTERFACE");
283
284 return init;
285}
Note: See TracBrowser for help on using the repository browser.