/*
* vmg - a versatile multigrid solver
* Copyright (C) 2012 Institute for Numerical Simulation, University of Bonn
*
* vmg is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* vmg is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see .
*/
/**
* @file mg.cpp
* @author Julian Iseringhausen
* @date Sat Jun 12 20:36:24 2010
*
* @brief A multigrid solver
*
* This file contains the implementation of the main class for
* a multigrid solver.
*
*/
#ifdef HAVE_CONFIG_H
#include
#endif
#ifdef DEBUG_MEASURE_TIME
#ifdef HAVE_MPI
#include
#ifdef HAVE_MARMOT
#include
#include
#endif
#endif
#endif
#include
#include
#include
#include
#include
#include
#include
#include "base/command_list.hpp"
#include "base/discretization.hpp"
#include "base/factory.hpp"
#include "base/interface.hpp"
#include "base/timer.hpp"
#include "comm/comm.hpp"
#include "grid/grid.hpp"
#include "grid/tempgrid.hpp"
#include "level/level_operator.hpp"
#include "smoother/smoother.hpp"
#include "solver/solver.hpp"
#include "mg.hpp"
using namespace VMG;
#define REGISTER_COMMAND(a) extern void Initialize##a();Initialize##a();
VMG::CommandFactory MG::command_factory;
static void VMGRegisterBuiltinCommands()
{
REGISTER_COMMAND(VMGCommandCheckConsistency);
REGISTER_COMMAND(VMGCommandCheckIterationCounter);
REGISTER_COMMAND(VMGCommandCheckRelativeResidual);
REGISTER_COMMAND(VMGCommandCheckResidual);
REGISTER_COMMAND(VMGCommandClearCoarseLevels);
REGISTER_COMMAND(VMGCommandClearGrid);
REGISTER_COMMAND(VMGCommandComputeResidualNorm);
REGISTER_COMMAND(VMGCommandCopyBoundary);
REGISTER_COMMAND(VMGCommandExecuteCycle);
REGISTER_COMMAND(VMGCommandExecuteCycleLoop);
REGISTER_COMMAND(VMGCommandExecuteFullCycle);
REGISTER_COMMAND(VMGCommandExecuteFullCycleLoop);
REGISTER_COMMAND(VMGCommandExportSolution);
REGISTER_COMMAND(VMGCommandForceDiscreteCompatibility);
REGISTER_COMMAND(VMGCommandImportRightHandSide);
REGISTER_COMMAND(VMGCommandInterpolateFMG);
REGISTER_COMMAND(VMGCommandInitializeIterationCounter);
REGISTER_COMMAND(VMGCommandInitializeResidualNorm);
REGISTER_COMMAND(VMGCommandNOP);
REGISTER_COMMAND(VMGCommandPrintAllSettings);
REGISTER_COMMAND(VMGCommandPrintDefect);
REGISTER_COMMAND(VMGCommandPrintGridStructure);
REGISTER_COMMAND(VMGCommandPrintGrid);
REGISTER_COMMAND(VMGCommandPrintResidualNorm);
REGISTER_COMMAND(VMGCommandPrintRunningTime);
REGISTER_COMMAND(VMGCommandProlongate);
REGISTER_COMMAND(VMGCommandRestrict);
REGISTER_COMMAND(VMGCommandSetAverageToZero);
REGISTER_COMMAND(VMGCommandSetCoarserDirichletValues);
REGISTER_COMMAND(VMGCommandSetLevel);
REGISTER_COMMAND(VMGCommandSmooth);
REGISTER_COMMAND(VMGCommandSolve);
}
MG::MG()
{
state = 0;
VMGRegisterBuiltinCommands();
}
MG::~MG()
{
MG::Destroy();
}
// Brings Multigrid back to starting state.
void MG::Destroy()
{
MG::Instance()->factories.clear();
MG::Instance()->state = 0;
Timer::Clear();
}
/*
* Post init communication class
*/
void MG::PostInit()
{
Multigrid* sol = new Multigrid(GetComm(), GetInterface());
sol->Register("SOL");
Multigrid* rhs = new Multigrid(GetComm(), GetInterface());
rhs->Register("RHS");
TempGrid* temp = new TempGrid();
temp->Register("TEMPGRID");
new ObjectStorage("GLOBAL_MAXLEVEL", sol->GlobalMaxLevel());
new ObjectStorage("MINLEVEL", sol->MinLevel());
new ObjectStorage("MAXLEVEL", sol->MaxLevel());
GetComm()->PostInit(*GetSol(), *GetRhs());
}
/**
* Solves a given system with a multigrid method
*
*/
void MG::Solve()
{
#ifdef DEBUG_MEASURE_TIME
#ifdef HAVE_MPI
GetComm()->Barrier();
#endif
Timer::Start("CompleteRunningTime");
#endif
CommandList* cl_init = MG::GetFactory().Get("COMMANDLIST_INIT")->Cast();
CommandList* cl_loop = MG::GetFactory().Get("COMMANDLIST_LOOP")->Cast();
CommandList* cl_finalize = MG::GetFactory().Get("COMMANDLIST_FINALIZE")->Cast();
cl_init->ExecuteList();
while (cl_loop->ExecuteList() == Continue);
cl_finalize->ExecuteList();
#ifdef DEBUG_MEASURE_TIME
#ifdef HAVE_MPI
GetComm()->Barrier();
#endif
Timer::Stop("CompleteRunningTime");
#ifdef DEBUG_MEASURE_TIME_OUTPUT
#ifdef HAVE_MPI
Timer::PrintGlobal();
#else
Timer::Print();
#endif
#endif
#endif
}
void MG::SetState(const int& state_)
{
MG::Instance()->state = state_;
}
VMG::Factory& MG::GetFactory()
{
std::map::iterator iter = MG::Instance()->factories.find(MG::Instance()->state);
if (iter == MG::Instance()->factories.end())
iter = MG::Instance()->factories.insert(std::make_pair(MG::Instance()->state, Factory())).first;
assert(iter != MG::Instance()->factories.end());
return iter->second;
}
VMG::CommandFactory& MG::GetCommands()
{
return MG::command_factory;
}
Comm* MG::GetComm()
{
return MG::GetFactory().Get("COMM")->Cast();
}
Discretization* MG::GetDiscretization()
{
return MG::GetFactory().Get("DISCRETIZATION")->Cast();
}
LevelOperator* MG::GetLevelOperator()
{
return MG::GetFactory().Get("LEVEL_OPERATOR")->Cast();
}
Multigrid* MG::GetRhs()
{
return MG::GetFactory().Get("RHS")->Cast();
}
Multigrid* MG::GetSol()
{
return MG::GetFactory().Get("SOL")->Cast();
}
VMG::Grid& MG::GetRhsMaxLevel()
{
return (*MG::GetRhs())(MG::GetRhs()->MaxLevel());
}
VMG::Grid& MG::GetSolMaxLevel()
{
return (*MG::GetSol())(MG::GetSol()->MaxLevel());
}
Smoother* MG::GetSmoother()
{
return MG::GetFactory().Get("SMOOTHER")->Cast();
}
Solver* MG::GetSolver()
{
return MG::GetFactory().Get("SOLVER")->Cast();
}
TempGrid* MG::GetTempGrid()
{
return MG::GetFactory().Get("TEMPGRID")->Cast();
}
Interface* MG::GetInterface()
{
return MG::GetFactory().Get("INTERFACE")->Cast();
}
static bool CheckObject(std::string id)
{
Object *obj = MG::GetFactory().Get(id);
#ifdef DEBUG_OUTPUT
if (obj == NULL)
printf("\nMultigrid: CLASS %s NOT INITIALIZED\n\n", id.c_str());
#endif
return obj != NULL;
}
bool MG::IsInitialized()
{
bool init = true;
init &= CheckObject("COMM");
init &= CheckObject("LEVEL_OPERATOR");
init &= CheckObject("RHS");
init &= CheckObject("SOL");
init &= CheckObject("SOLVER");
init &= CheckObject("SMOOTHER");
init &= CheckObject("DISCRETIZATION");
init &= CheckObject("MAX_ITERATION");
init &= CheckObject("PRECISION");
init &= CheckObject("PRESMOOTHSTEPS");
init &= CheckObject("POSTSMOOTHSTEPS");
init &= CheckObject("COMMANDLIST_INIT");
init &= CheckObject("COMMANDLIST_LOOP");
init &= CheckObject("COMMANDLIST_FINALIZE");
init &= CheckObject("MINLEVEL");
init &= CheckObject("MAXLEVEL");
init &= CheckObject("GLOBAL_MAXLEVEL");
init &= CheckObject("INTERFACE");
return init;
}