Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
Nektar::MultiRegions::GlobalLinSysIterativeFull Class Reference

A global linear system. More...

#include <GlobalLinSysIterativeFull.h>

Inheritance diagram for Nektar::MultiRegions::GlobalLinSysIterativeFull:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::GlobalLinSysIterativeFull:
Collaboration graph
[legend]

Public Member Functions

 GlobalLinSysIterativeFull (const GlobalLinSysKey &pLinSysKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve.
virtual ~GlobalLinSysIterativeFull ()
- Public Member Functions inherited from Nektar::MultiRegions::GlobalLinSysIterative
 GlobalLinSysIterative (const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve.
virtual ~GlobalLinSysIterative ()
- Public Member Functions inherited from Nektar::MultiRegions::GlobalLinSys
 GlobalLinSys (const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve.
virtual ~GlobalLinSys ()
const GlobalLinSysKeyGetKey (void) const
 Returns the key associated with the system.
const boost::weak_ptr< ExpList > & GetLocMat (void) const
void InitObject ()
void Initialise (const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
void Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve the linear system for given input and output vectors using a specified local to global map.
boost::shared_ptr< GlobalLinSysGetSharedThisPtr ()
 Returns a shared pointer to the current object.
int GetNumBlocks ()
DNekScalMatSharedPtr GetBlock (unsigned int n)
DNekScalBlkMatSharedPtr GetStaticCondBlock (unsigned int n)
void DropStaticCondBlock (unsigned int n)
void SolveLinearSystem (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir=0)
 Solve the linear system for given input and output vectors.

Static Public Member Functions

static GlobalLinSysSharedPtr create (const GlobalLinSysKey &pLinSysKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Creates an instance of this class.

Static Public Attributes

static std::string className
 Name of class.

Private Member Functions

virtual void v_Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve the linear system for given input and output vectors using a specified local to global map.
virtual void v_DoMatrixMultiply (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
virtual void v_UniqueMap ()

Private Attributes

boost::shared_ptr< AssemblyMapm_locToGloMap

Additional Inherited Members

- Protected Member Functions inherited from Nektar::MultiRegions::GlobalLinSysIterative
void DoAconjugateProjection (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
 A-conjugate projection technique.
void DoConjugateGradient (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)
 Actual iterative solve.
void Set_Rhs_Magnitude (const NekVector< NekDouble > &pIn)
virtual void v_UniqueMap ()=0
- Protected Attributes inherited from Nektar::MultiRegions::GlobalLinSysIterative
Array< OneD, int > m_map
 Global to universal unique map.
NekDouble m_tolerance
 Tolerance of iterative solver.
NekDouble m_rhs_magnitude
 dot product of rhs to normalise stopping criterion
PreconditionerSharedPtr m_precon
MultiRegions::PreconditionerType m_precontype
int m_totalIterations
bool m_useProjection
 Whether to apply projection technique.
bool m_root
 Provide verbose output and root if parallel.
bool m_verbose
boost::circular_buffer< Array
< OneD, NekDouble > > 
m_prevLinSol
 Storage for solutions to previous linear problems.
int m_numPrevSols
 Total counter of previous solutions.

Detailed Description

A global linear system.

Definition at line 48 of file GlobalLinSysIterativeFull.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::GlobalLinSysIterativeFull::GlobalLinSysIterativeFull ( const GlobalLinSysKey pKey,
const boost::weak_ptr< ExpList > &  pExp,
const boost::shared_ptr< AssemblyMap > &  pLocToGloMap 
)

Constructor for full direct matrix solve.

Constructor for full direct matrix solve.

Parameters
pKeyKey specifying matrix to solve.
pExpShared pointer to expansion list for applying matrix evaluations.
pLocToGloMapLocal to global mapping.

Definition at line 69 of file GlobalLinSysIterativeFull.cpp.

References ASSERTL1, Nektar::MultiRegions::eIterativeFull, Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType(), and Nektar::MultiRegions::GlobalLinSys::m_linSysKey.

: GlobalLinSys (pKey, pExp, pLocToGloMap),
GlobalLinSysIterative(pKey, pExp, pLocToGloMap)
{
"This routine should only be used when using an Iterative "
"conjugate gradient matrix solve.");
}
Nektar::MultiRegions::GlobalLinSysIterativeFull::~GlobalLinSysIterativeFull ( )
virtual

Definition at line 85 of file GlobalLinSysIterativeFull.cpp.

{
}

Member Function Documentation

static GlobalLinSysSharedPtr Nektar::MultiRegions::GlobalLinSysIterativeFull::create ( const GlobalLinSysKey pLinSysKey,
const boost::weak_ptr< ExpList > &  pExpList,
const boost::shared_ptr< AssemblyMap > &  pLocToGloMap 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file GlobalLinSysIterativeFull.h.

{
::AllocateSharedPtr(pLinSysKey, pExpList, pLocToGloMap);
}
void Nektar::MultiRegions::GlobalLinSysIterativeFull::v_DoMatrixMultiply ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Definition at line 183 of file GlobalLinSysIterativeFull.cpp.

References ASSERTL0, Nektar::MultiRegions::eGlobal, Nektar::eWrapper, Nektar::MultiRegions::GlobalLinSys::m_expList, Nektar::MultiRegions::GlobalLinSys::m_linSysKey, m_locToGloMap, Vmath::Vadd(), and Vmath::Zero().

{
boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
// Perform matrix-vector operation A*d_i
expList->GeneralMatrixOp(m_linSysKey,
pInput, pOutput, eGlobal);
// retrieve robin boundary condition information and apply robin
// boundary conditions to the solution.
const std::map<int, RobinBCInfoSharedPtr> vRobinBCInfo
= expList->GetRobinBCInfo();
if(vRobinBCInfo.size() > 0)
{
ASSERTL0(false,
"Robin boundaries not set up in IterativeFull solver.");
int nGlobal = m_locToGloMap->GetNumGlobalCoeffs();
int nLocal = m_locToGloMap->GetNumLocalCoeffs();
int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
int nNonDir = nGlobal - nDir;
Array<OneD, NekDouble> robin_A(nGlobal, 0.0);
Array<OneD, NekDouble> robin_l(nLocal, 0.0);
Array<OneD, NekDouble> tmp;
NekVector<NekDouble> robin(nNonDir,
tmp = robin_A + nDir, eWrapper);
// Operation: p_A = A * d_A
// First map d_A to local solution
m_locToGloMap->GlobalToLocal(pInput, robin_l);
// Iterate over all the elements computing Robin BCs where
// necessary
for (int n = 0; n < expList->GetNumElmts(); ++n)
{
int nel = expList->GetOffset_Elmt_Id(n);
int offset = expList->GetCoeff_Offset(n);
int ncoeffs = expList->GetExp(nel)->GetNcoeffs();
if(vRobinBCInfo.count(nel) != 0) // add robin mass matrix
{
Array<OneD, NekDouble> tmp;
StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(nel);
// add local matrix contribution
for(rBC = vRobinBCInfo.find(nel)->second;rBC; rBC = rBC->next)
{
vExp->AddRobinEdgeContribution(rBC->m_robinID,rBC->m_robinPrimitiveCoeffs, tmp = robin_l + offset);
}
}
else
{
Vmath::Zero(ncoeffs, &robin_l[offset], 1);
}
}
// Map local Robin contribution back to global coefficients
m_locToGloMap->LocalToGlobal(robin_l, robin_A);
// Add them to the output of the GeneralMatrixOp
Vmath::Vadd(nGlobal, pOutput, 1, robin_A, 1, pOutput, 1);
}
}
void Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve ( const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const AssemblyMapSharedPtr pLocToGloMap,
const Array< OneD, const NekDouble > &  pDirForcing = NullNekDouble1DArray 
)
privatevirtual

Solve the linear system for given input and output vectors using a specified local to global map.

Solve a global linear system with Dirichlet forcing using a conjugate gradient method. This routine performs handling of the Dirichlet forcing terms and wraps the underlying iterative solver used for the remaining degrees of freedom.

Consider solving for $x$, the matrix system $Ax=b$, where $b$ is known. To enforce the Dirichlet terms we instead solve

\[A(x-x_0) = b - Ax_0 \]

where $x_0$ is the Dirichlet forcing.

Parameters
pInputRHS of linear system, $b$.
pOutputOn input, values of dirichlet degrees of freedom with initial guess on other values. On output, the solution $x$.
pLocToGloMapLocal to global mapping.
pDirForcingPrecalculated Dirichlet forcing.

Definition at line 109 of file GlobalLinSysIterativeFull.cpp.

References ASSERTL0, Nektar::MultiRegions::eGlobal, Nektar::MultiRegions::GlobalLinSys::m_expList, Nektar::MultiRegions::GlobalLinSys::m_linSysKey, m_locToGloMap, Nektar::LibUtilities::ReduceSum, Nektar::MultiRegions::GlobalLinSys::SolveLinearSystem(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vsub().

{
boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
bool vCG;
if ((m_locToGloMap = boost::dynamic_pointer_cast<AssemblyMapCG>(
pLocToGloMap)))
{
vCG = true;
}
else if ((m_locToGloMap = boost::dynamic_pointer_cast<
AssemblyMapDG>(pLocToGloMap)))
{
vCG = false;
}
bool dirForcCalculated = (bool) pDirForcing.num_elements();
int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
int nDirTotal = nDirDofs;
expList->GetComm()->AllReduce(nDirTotal, LibUtilities::ReduceSum);
Array<OneD, NekDouble> tmp(nGlobDofs), tmp2;
if(nDirTotal)
{
// calculate the Dirichlet forcing
if(dirForcCalculated)
{
Vmath::Vsub(nGlobDofs, pInput.get(), 1,
pDirForcing.get(), 1,
tmp.get(), 1);
}
else
{
// Calculate the dirichlet forcing B_b (== X_b) and
// substract it from the rhs
expList->GeneralMatrixOp(
m_linSysKey, pOutput, tmp, eGlobal);
Vmath::Vsub(nGlobDofs, pInput.get(), 1,
tmp.get(), 1,
tmp.get(), 1);
}
if (vCG)
{
Array<OneD, NekDouble> out(nGlobDofs,0.0);
// solve for perturbation from intiial guess in pOutput
nGlobDofs, tmp, out, pLocToGloMap, nDirDofs);
Vmath::Vadd(nGlobDofs-nDirDofs, &out [nDirDofs], 1,
&pOutput[nDirDofs], 1, &pOutput[nDirDofs], 1);
}
else
{
ASSERTL0(false, "Need DG solve if using Dir BCs");
}
}
else
{
Vmath::Vcopy(nGlobDofs, pInput, 1, tmp, 1);
SolveLinearSystem(nGlobDofs, tmp, pOutput, pLocToGloMap);
}
}
void Nektar::MultiRegions::GlobalLinSysIterativeFull::v_UniqueMap ( )
privatevirtual

Definition at line 251 of file GlobalLinSysIterativeFull.cpp.

References m_locToGloMap, and Nektar::MultiRegions::GlobalLinSysIterative::m_map.

{
m_map = m_locToGloMap->GetGlobalToUniversalMapUnique();
}

Member Data Documentation

string Nektar::MultiRegions::GlobalLinSysIterativeFull::className
static
Initial value:
"IterativeFull",
"Iterative solver for full matrix system.")

Name of class.

Registers the class with the Factory.

Definition at line 63 of file GlobalLinSysIterativeFull.h.

boost::shared_ptr<AssemblyMap> Nektar::MultiRegions::GlobalLinSysIterativeFull::m_locToGloMap
private

Definition at line 76 of file GlobalLinSysIterativeFull.h.

Referenced by v_DoMatrixMultiply(), v_Solve(), and v_UniqueMap().