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

A global linear system. More...

#include <GlobalLinSys.h>

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

Public Member Functions

 GlobalLinSys (const GlobalLinSysKey &pKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve. More...
 
virtual ~GlobalLinSys ()
 
const GlobalLinSysKeyGetKey (void) const
 Returns the key associated with the system. More...
 
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. More...
 
boost::shared_ptr< GlobalLinSysGetSharedThisPtr ()
 Returns a shared pointer to the current object. More...
 
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. More...
 

Protected Member Functions

virtual int v_GetNumBlocks ()
 Get the number of blocks in this system. More...
 
virtual DNekScalMatSharedPtr v_GetBlock (unsigned int n)
 Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock (unsigned int n)
 Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual void v_DropStaticCondBlock (unsigned int n)
 Releases the static condensation block matrices from NekManager of n-th expansion using the matrix key provided by the m_linSysKey. More...
 
PreconditionerSharedPtr CreatePrecon (AssemblyMapSharedPtr asmMap)
 Create a preconditioner object from the parameters defined in the supplied assembly map. More...
 

Protected Attributes

const GlobalLinSysKey m_linSysKey
 Key associated with this linear system. More...
 
const boost::weak_ptr< ExpListm_expList
 Local Matrix System. More...
 
const std::map< int,
RobinBCInfoSharedPtr
m_robinBCInfo
 Robin boundary info. More...
 

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)=0
 Solve a linear system based on mapping. More...
 
virtual void v_SolveLinearSystem (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)=0
 Solve a basic matrix system. More...
 
virtual void v_InitObject ()
 
virtual void v_Initialise (const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 

Static Private Attributes

static std::string lookupIds []
 
static std::string def
 

Detailed Description

A global linear system.

Consider the linear system $\boldsymbol{M\hat{u}}_g=\boldsymbol{\hat{f}}$. Distinguishing between the boundary and interior components of $\boldsymbol{\hat{u}}_g$ and $\boldsymbol{\hat{f}}$ using $\boldsymbol{\hat{u}}_b$, $\boldsymbol{\hat{u}}_i$ and $\boldsymbol{\hat{f}}_b$, $\boldsymbol{\hat{f}}_i$ respectively, this system can be split into its constituent parts as

\[\left[\begin{array}{cc} \boldsymbol{M}_b&\boldsymbol{M}_{c1}\\ \boldsymbol{M}_{c2}&\boldsymbol{M}_i\\ \end{array}\right] \left[\begin{array}{c} \boldsymbol{\hat{u}_b}\\ \boldsymbol{\hat{u}_i}\\ \end{array}\right]= \left[\begin{array}{c} \boldsymbol{\hat{f}_b}\\ \boldsymbol{\hat{f}_i}\\ \end{array}\right]\]

where $\boldsymbol{M}_b$ represents the components of $\boldsymbol{M}$ resulting from boundary-boundary mode interactions, $\boldsymbol{M}_{c1}$ and $\boldsymbol{M}_{c2}$ represent the components resulting from coupling between the boundary-interior modes, and $\boldsymbol{M}_i$ represents the components of $\boldsymbol{M}$ resulting from interior-interior mode interactions.

The solution of the linear system can now be determined in two steps:

\begin{eqnarray*} \mathrm{step 1:}&\quad&(\boldsymbol{M}_b-\boldsymbol{M}_{c1} \boldsymbol{M}_i^{-1}\boldsymbol{M}_{c2}) \boldsymbol{\hat{u}_b} = \boldsymbol{\hat{f}}_b - \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1} \boldsymbol{\hat{f}}_i,\nonumber \\ \mathrm{step 2:}&\quad&\boldsymbol{\hat{u}_i}=\boldsymbol{M}_i^{-1} \left( \boldsymbol{\hat{f}}_i - \boldsymbol{M}_{c2}\boldsymbol{\hat{u}_b} \right). \nonumber \\ \end{eqnarray*}

As the inverse of $\boldsymbol{M}_i^{-1}$ is

\[ \boldsymbol{M}_i^{-1} = \left [\underline{\boldsymbol{M}^e_i} \right ]^{-1} = \underline{[\boldsymbol{M}^e_i]}^{-1} \]

and the following operations can be evaluated as,

\begin{eqnarray*} \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}\boldsymbol{\hat{f}}_i & =& \boldsymbol{\mathcal{A}}_b^T \underline{\boldsymbol{M}^e_{c1}} \underline{[\boldsymbol{M}^e_i]}^{-1} \boldsymbol{\hat{f}}_i \\ \boldsymbol{M}_{c2} \boldsymbol{\hat{u}_b} &=& \underline{\boldsymbol{M}^e_{c2}} \boldsymbol{\mathcal{A}}_b \boldsymbol{\hat{u}_b}.\end{eqnarray*}

where $\boldsymbol{\mathcal{A}}_b $ is the permutation matrix which scatters from global to local degrees of freedom, only the following four matrices should be constructed:

The first three matrices are just a concatenation of the corresponding local matrices and they can be created as such. They also allow for an elemental evaluation of the operations concerned.

The global Schur complement however should be assembled from the concatenation of the local elemental Schur complements, that is,

\[ \boldsymbol{M}_{\mathrm{Schur}}=\boldsymbol{M}_b - \boldsymbol{M}_{c1} \boldsymbol{M}_i^{-1} \boldsymbol{M}_{c2} = \boldsymbol{\mathcal{A}}_b^T \left [\underline{\boldsymbol{M}^e_b - \boldsymbol{M}^e_{c1} [\boldsymbol{M}^e_i]^{-1} (\boldsymbol{M}^e_{c2})} \right ] \boldsymbol{\mathcal{A}}_b \]

and it is the only matrix operation that need to be evaluated on a global level when using static condensation. However, due to the size and sparsity of the matrix $\boldsymbol{\mathcal{A}}_b$, it is more efficient to assemble the global Schur matrix using the mapping array bmap $[e][i]$ contained in the input argument locToGloMap. The global Schur complement is then constructed as:

\[\boldsymbol{M}_{\mathrm{Schur}}\left[\mathrm{\a bmap}[e][i]\right] \left[\mathrm{\a bmap}[e][j]\right]=\mathrm{\a bsign}[e][i]\cdot \mathrm{\a bsign}[e][j] \cdot\boldsymbol{M}^e_{\mathrm{Schur}}[i][j]\]

All four matrices are stored in the GlobalLinSys returned by this function.

Definition at line 74 of file GlobalLinSys.h.

Constructor & Destructor Documentation

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.

Given a block matrix, construct a global matrix system according to a local to global mapping. #m_linSys is constructed by AssembleFullMatrix().

Parameters
pkeyAssociated linear system key.
locToGloMapLocal to global mapping.

Definition at line 189 of file GlobalLinSys.cpp.

192  :
193  m_linSysKey(pKey),
194  m_expList(pExpList),
195  m_robinBCInfo(m_expList.lock()->GetRobinBCInfo())
196  {
197  }
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:127
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:131
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
virtual Nektar::MultiRegions::GlobalLinSys::~GlobalLinSys ( )
inlinevirtual

Definition at line 84 of file GlobalLinSys.h.

84 {}

Member Function Documentation

PreconditionerSharedPtr Nektar::MultiRegions::GlobalLinSys::CreatePrecon ( AssemblyMapSharedPtr  asmMap)
protected

Create a preconditioner object from the parameters defined in the supplied assembly map.

Parameters
asmMapAssembly map used to construct the global system.

Definition at line 217 of file GlobalLinSys.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::MultiRegions::GetPreconFactory(), GetSharedThisPtr(), and Nektar::MultiRegions::PreconditionerTypeMap.

Referenced by Nektar::MultiRegions::GlobalLinSysIterative::DoConjugateGradient(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_PreSolve(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_PreSolve(), and Nektar::MultiRegions::GlobalLinSysPETSc::v_SolveLinearSystem().

218  {
219  PreconditionerType pType = asmMap->GetPreconType();
220  std::string PreconType = MultiRegions::PreconditionerTypeMap[pType];
222  PreconType, GetSharedThisPtr(), asmMap);
223  }
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
PreconFactory & GetPreconFactory()
boost::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:107
const char *const PreconditionerTypeMap[]
void Nektar::MultiRegions::GlobalLinSys::DropStaticCondBlock ( unsigned int  n)
inline

Definition at line 230 of file GlobalLinSys.h.

References v_DropStaticCondBlock().

231  {
232  return v_DropStaticCondBlock(n);
233  }
virtual void v_DropStaticCondBlock(unsigned int n)
Releases the static condensation block matrices from NekManager of n-th expansion using the matrix ke...
DNekScalMatSharedPtr Nektar::MultiRegions::GlobalLinSys::GetBlock ( unsigned int  n)
inline

Definition at line 220 of file GlobalLinSys.h.

References v_GetBlock().

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::AssembleFullMatrix(), Nektar::MultiRegions::GlobalLinSysXxtFull::AssembleMatrixArrays(), and Nektar::MultiRegions::GlobalLinSysPETScFull::GlobalLinSysPETScFull().

221  {
222  return v_GetBlock(n);
223  }
virtual DNekScalMatSharedPtr v_GetBlock(unsigned int n)
Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey...
const GlobalLinSysKey & Nektar::MultiRegions::GlobalLinSys::GetKey ( void  ) const
inline

Returns the key associated with the system.

Definition at line 169 of file GlobalLinSys.h.

References m_linSysKey.

170  {
171  return m_linSysKey;
172  }
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:127
const boost::weak_ptr< ExpList > & Nektar::MultiRegions::GlobalLinSys::GetLocMat ( void  ) const
inline

Definition at line 177 of file GlobalLinSys.h.

References m_expList.

178  {
179  return m_expList;
180  }
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
int Nektar::MultiRegions::GlobalLinSys::GetNumBlocks ( )
inline

Definition at line 235 of file GlobalLinSys.h.

References v_GetNumBlocks().

236  {
237  return v_GetNumBlocks();
238  }
virtual int v_GetNumBlocks()
Get the number of blocks in this system.
boost::shared_ptr<GlobalLinSys> Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr ( )
inline

Returns a shared pointer to the current object.

Definition at line 107 of file GlobalLinSys.h.

Referenced by CreatePrecon().

108  {
109  return shared_from_this();
110  }
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSys::GetStaticCondBlock ( unsigned int  n)
inline

Definition at line 225 of file GlobalLinSys.h.

References v_GetStaticCondBlock().

226  {
227  return v_GetStaticCondBlock(n);
228  }
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock(unsigned int n)
Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided ...
void Nektar::MultiRegions::GlobalLinSys::Initialise ( const boost::shared_ptr< AssemblyMap > &  pLocToGloMap)
inline
void Nektar::MultiRegions::GlobalLinSys::InitObject ( )
inline

Definition at line 209 of file GlobalLinSys.h.

References v_InitObject().

210  {
211  v_InitObject();
212  }
void Nektar::MultiRegions::GlobalLinSys::Solve ( const Array< OneD, const NekDouble > &  in,
Array< OneD, NekDouble > &  out,
const AssemblyMapSharedPtr locToGloMap,
const Array< OneD, const NekDouble > &  dirForcing = NullNekDouble1DArray 
)
inline

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

Definition at line 186 of file GlobalLinSys.h.

References v_Solve().

191  {
192  v_Solve(in,out,locToGloMap,dirForcing);
193  }
virtual void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)=0
Solve a linear system based on mapping.
void Nektar::MultiRegions::GlobalLinSys::SolveLinearSystem ( const int  pNumRows,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const AssemblyMapSharedPtr locToGloMap,
const int  pNumDir = 0 
)
inline

Solve the linear system for given input and output vectors.

Definition at line 199 of file GlobalLinSys.h.

References v_SolveLinearSystem().

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysXxtFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysPETScFull::v_Solve(), and Nektar::MultiRegions::GlobalLinSysStaticCond::v_Solve().

205  {
206  v_SolveLinearSystem(pNumRows, pInput, pOutput, locToGloMap, pNumDir);
207  }
virtual void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)=0
Solve a basic matrix system.
void Nektar::MultiRegions::GlobalLinSys::v_DropStaticCondBlock ( unsigned int  n)
protectedvirtual

Releases the static condensation block matrices from NekManager of n-th expansion using the matrix key provided by the m_linSysKey.

Parameters
nNumber of the expansion

Definition at line 394 of file GlobalLinSys.cpp.

References Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_expList, and m_linSysKey.

Referenced by DropStaticCondBlock(), and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::PrepareLocalSchurComplement().

395  {
396  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
397 
398  StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(n);
399 
400  // need to be initialised with zero size for non variable
401  // coefficient case
402  StdRegions::VarCoeffMap vVarCoeffMap;
403 
404  // retrieve variable coefficients
405  if(m_linSysKey.GetNVarCoeffs() > 0)
406  {
407  StdRegions::VarCoeffMap::const_iterator x;
408  int cnt = expList->GetPhys_Offset(n);
409  for (x = m_linSysKey.GetVarCoeffs().begin();
410  x != m_linSysKey.GetVarCoeffs().end (); ++x)
411  {
412  vVarCoeffMap[x->first] = x->second + cnt;
413  }
414  }
415 
416  LocalRegions::MatrixKey matkey(m_linSysKey.GetMatrixType(),
417  vExp->DetShapeType(),
418  *vExp,
420  vVarCoeffMap);
421 
422  vExp->DropLocStaticCondMatrix(matkey);
423  }
const StdRegions::VarCoeffMap & GetVarCoeffs() const
const StdRegions::ConstFactorMap & GetConstFactors() const
Returns all the constants.
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:127
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
DNekScalMatSharedPtr Nektar::MultiRegions::GlobalLinSys::v_GetBlock ( unsigned int  n)
protectedvirtual

Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey.

Parameters
nNumber of the expansion.
Returns
Block matrix for the specified expansion.

Definition at line 244 of file GlobalLinSys.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_expList, m_linSysKey, and m_robinBCInfo.

Referenced by GetBlock(), and Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel().

245  {
246  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
247  int cnt = 0;
248  DNekScalMatSharedPtr loc_mat;
249 
251  boost::dynamic_pointer_cast<LocalRegions::Expansion>(
252  expList->GetExp(n));
253 
254  // need to be initialised with zero size for non variable
255  // coefficient case
256  StdRegions::VarCoeffMap vVarCoeffMap;
257 
258  // retrieve variable coefficients
259  if(m_linSysKey.GetNVarCoeffs() > 0)
260  {
261  StdRegions::VarCoeffMap::const_iterator x;
262  cnt = expList->GetPhys_Offset(n);
263 
264  for (x = m_linSysKey.GetVarCoeffs().begin();
265  x != m_linSysKey.GetVarCoeffs().end(); ++x)
266  {
267  vVarCoeffMap[x->first] = x->second + cnt;
268  }
269  }
270 
271  LocalRegions::MatrixKey matkey(m_linSysKey.GetMatrixType(),
272  vExp->DetShapeType(),
273  *vExp, m_linSysKey.GetConstFactors(),
274  vVarCoeffMap);
275  loc_mat = vExp->GetLocMatrix(matkey);
276 
277  // apply robin boundary conditions to the matrix.
278  if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
279  {
281 
282  // declare local matrix from scaled matrix.
283  int rows = loc_mat->GetRows();
284  int cols = loc_mat->GetColumns();
285  const NekDouble *dat = loc_mat->GetRawPtr();
287  AllocateSharedPtr(rows,cols,dat);
288  Blas::Dscal(rows*cols,loc_mat->Scale(),new_mat->GetRawPtr(),1);
289 
290  // add local matrix contribution
291  for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
292  {
293  vExp->AddRobinMassMatrix(
294  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
295  }
296 
297  // redeclare loc_mat to point to new_mat plus the scalar.
299  1.0, new_mat);
300  }
301 
302  // finally return the matrix.
303  return loc_mat;
304  }
const StdRegions::VarCoeffMap & GetVarCoeffs() const
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
const StdRegions::ConstFactorMap & GetConstFactors() const
Returns all the constants.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
double NekDouble
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:127
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:131
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
int Nektar::MultiRegions::GlobalLinSys::v_GetNumBlocks ( )
protectedvirtual

Get the number of blocks in this system.

At the top level this corresponds to the number of elements in the expansion list.

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 232 of file GlobalLinSys.cpp.

References m_expList.

Referenced by GetNumBlocks().

233  {
234  return m_expList.lock()->GetExpSize();
235  }
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSys::v_GetStaticCondBlock ( unsigned int  n)
protectedvirtual

Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided by the m_linSysKey.

Parameters
nNumber of the expansion
Returns
2x2 Block matrix holding the static condensation matrices for the n-th expansion.

Reimplemented in Nektar::MultiRegions::GlobalLinSysIterativeStaticCond, and Nektar::MultiRegions::GlobalLinSysPETScStaticCond.

Definition at line 314 of file GlobalLinSys.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::MultiRegions::GlobalMatrixKey::GetConstFactors(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalMatrixKey::GetNVarCoeffs(), Nektar::MultiRegions::GlobalMatrixKey::GetVarCoeffs(), m_expList, m_linSysKey, and m_robinBCInfo.

Referenced by GetStaticCondBlock(), and Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel().

316  {
317  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
318  int cnt = 0;
319  DNekScalBlkMatSharedPtr loc_mat;
320  DNekScalMatSharedPtr tmp_mat;
321 
322  StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(n);
323 
324  // need to be initialised with zero size for non variable
325  // coefficient case
326  StdRegions::VarCoeffMap vVarCoeffMap;
327 
328  // retrieve variable coefficients
329  if(m_linSysKey.GetNVarCoeffs() > 0)
330  {
331  StdRegions::VarCoeffMap::const_iterator x;
332  cnt = expList->GetPhys_Offset(n);
333  for (x = m_linSysKey.GetVarCoeffs().begin();
334  x != m_linSysKey.GetVarCoeffs().end (); ++x)
335  {
336  vVarCoeffMap[x->first] = x->second + cnt;
337  }
338  }
339 
340  LocalRegions::MatrixKey matkey(m_linSysKey.GetMatrixType(),
341  vExp->DetShapeType(),
342  *vExp,
344  vVarCoeffMap);
345 
346  loc_mat = vExp->GetLocStaticCondMatrix(matkey);
347 
348  if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
349  {
351 
352  tmp_mat = loc_mat->GetBlock(0,0);
353 
354  // declare local matrix from scaled matrix.
355  int rows = tmp_mat->GetRows();
356  int cols = tmp_mat->GetColumns();
357  const NekDouble *dat = tmp_mat->GetRawPtr();
359  AllocateSharedPtr(rows, cols, dat);
360  Blas::Dscal(rows*cols,tmp_mat->Scale(),new_mat->GetRawPtr(),1);
361 
362  // add local matrix contribution
363  for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
364  {
365  vExp->AddRobinMassMatrix(
366  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
367  }
368 
369  // redeclare loc_mat to point to new_mat plus the scalar.
371  1.0, new_mat);
372  DNekScalBlkMatSharedPtr new_loc_mat;
373  unsigned int exp_size[] = {tmp_mat->GetRows(), loc_mat->GetBlock(1,1)->GetRows()};
374  unsigned int nblks = 2;
375  new_loc_mat = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size);
376 
377 
378  new_loc_mat->SetBlock(0,0,tmp_mat);
379  new_loc_mat->SetBlock(0,1,loc_mat->GetBlock(0,1));
380  new_loc_mat->SetBlock(1,0,loc_mat->GetBlock(1,0));
381  new_loc_mat->SetBlock(1,1,loc_mat->GetBlock(1,1));
382  loc_mat = new_loc_mat;
383  }
384 
385  return loc_mat;
386  }
const StdRegions::VarCoeffMap & GetVarCoeffs() const
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
const StdRegions::ConstFactorMap & GetConstFactors() const
Returns all the constants.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
double NekDouble
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:127
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:131
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129
void Nektar::MultiRegions::GlobalLinSys::v_Initialise ( const boost::shared_ptr< AssemblyMap > &  pLocToGloMap)
privatevirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 430 of file GlobalLinSys.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by Initialise().

432  {
433  NEKERROR(ErrorUtil::efatal, "Method does not exist" );
434  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:185
void Nektar::MultiRegions::GlobalLinSys::v_InitObject ( )
privatevirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysIterativeStaticCond, Nektar::MultiRegions::GlobalLinSysStaticCond, and Nektar::MultiRegions::GlobalLinSysPETScStaticCond.

Definition at line 425 of file GlobalLinSys.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by InitObject().

426  {
427  NEKERROR(ErrorUtil::efatal, "Method does not exist" );
428  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:185
virtual void Nektar::MultiRegions::GlobalLinSys::v_Solve ( const Array< OneD, const NekDouble > &  in,
Array< OneD, NekDouble > &  out,
const AssemblyMapSharedPtr locToGloMap,
const Array< OneD, const NekDouble > &  dirForcing = NullNekDouble1DArray 
)
privatepure virtual
virtual void Nektar::MultiRegions::GlobalLinSys::v_SolveLinearSystem ( const int  pNumRows,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const AssemblyMapSharedPtr locToGloMap,
const int  pNumDir 
)
privatepure virtual

Member Data Documentation

std::string Nektar::MultiRegions::GlobalLinSys::def
staticprivate
Initial value:
"DirectMultiLevelStaticCond")

Definition at line 162 of file GlobalLinSys.h.

std::string Nektar::MultiRegions::GlobalLinSys::lookupIds
staticprivate

Definition at line 161 of file GlobalLinSys.h.

const boost::weak_ptr<ExpList> Nektar::MultiRegions::GlobalLinSys::m_expList
protected

Local Matrix System.

Definition at line 129 of file GlobalLinSys.h.

Referenced by Nektar::MultiRegions::GlobalLinSysDirectFull::AssembleFullMatrix(), Nektar::MultiRegions::GlobalLinSysXxtFull::AssembleMatrixArrays(), Nektar::MultiRegions::GlobalLinSysIterative::CalculateAnorm(), Nektar::MultiRegions::GlobalLinSysPETSc::CalculateReordering(), Nektar::MultiRegions::GlobalLinSysStaticCond::ConstructNextLevelCondensedSystem(), Nektar::MultiRegions::GlobalLinSysIterative::DoAconjugateProjection(), Nektar::MultiRegions::GlobalLinSysIterative::DoConjugateGradient(), GetLocMat(), Nektar::MultiRegions::GlobalLinSysIterative::GlobalLinSysIterative(), Nektar::MultiRegions::GlobalLinSysPETScFull::GlobalLinSysPETScFull(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::PrepareLocalSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterative::Set_Rhs_Magnitude(), Nektar::MultiRegions::GlobalLinSysStaticCond::SetupTopLevel(), Nektar::MultiRegions::GlobalLinSysIterative::UpdateKnownSolutions(), Nektar::MultiRegions::GlobalLinSysXxtStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_AssembleSchurComplement(), Nektar::MultiRegions::GlobalLinSysPETScFull::v_DoMatrixMultiply(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_DoMatrixMultiply(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_DoMatrixMultiply(), v_DropStaticCondBlock(), v_GetBlock(), v_GetNumBlocks(), v_GetStaticCondBlock(), Nektar::MultiRegions::GlobalLinSysPETScStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject(), Nektar::MultiRegions::GlobalLinSysDirectFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysXxtFull::v_Solve(), Nektar::MultiRegions::GlobalLinSysIterativeFull::v_Solve(), and Nektar::MultiRegions::GlobalLinSysPETScFull::v_Solve().

const GlobalLinSysKey Nektar::MultiRegions::GlobalLinSys::m_linSysKey
protected
const std::map<int, RobinBCInfoSharedPtr> Nektar::MultiRegions::GlobalLinSys::m_robinBCInfo
protected

Robin boundary info.

Definition at line 131 of file GlobalLinSys.h.

Referenced by v_GetBlock(), and v_GetStaticCondBlock().