Nektar++
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...
 

Protected Attributes

const GlobalLinSysKey m_linSysKey
 Key associated with this linear system. More...
 
const boost::weak_ptr< ExpListm_expList
 Local Matrix System. More...
 
const map< int, RobinBCInfoSharedPtrm_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 70 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 188 of file GlobalLinSys.cpp.

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

Definition at line 80 of file GlobalLinSys.h.

80 {}

Member Function Documentation

void Nektar::MultiRegions::GlobalLinSys::DropStaticCondBlock ( unsigned int  n)
inline

Definition at line 224 of file GlobalLinSys.h.

References v_DropStaticCondBlock().

225  {
226  return v_DropStaticCondBlock(n);
227  }
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 214 of file GlobalLinSys.h.

References v_GetBlock().

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

215  {
216  return v_GetBlock(n);
217  }
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 163 of file GlobalLinSys.h.

References m_linSysKey.

164  {
165  return m_linSysKey;
166  }
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:123
const boost::weak_ptr< ExpList > & Nektar::MultiRegions::GlobalLinSys::GetLocMat ( void  ) const
inline

Definition at line 171 of file GlobalLinSys.h.

References m_expList.

172  {
173  return m_expList;
174  }
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
int Nektar::MultiRegions::GlobalLinSys::GetNumBlocks ( )
inline

Definition at line 229 of file GlobalLinSys.h.

References v_GetNumBlocks().

230  {
231  return v_GetNumBlocks();
232  }
virtual int v_GetNumBlocks()
Get the number of blocks in this system.
boost::shared_ptr<GlobalLinSys> Nektar::MultiRegions::GlobalLinSys::GetSharedThisPtr ( )
inline
DNekScalBlkMatSharedPtr Nektar::MultiRegions::GlobalLinSys::GetStaticCondBlock ( unsigned int  n)
inline

Definition at line 219 of file GlobalLinSys.h.

References v_GetStaticCondBlock().

220  {
221  return v_GetStaticCondBlock(n);
222  }
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

Definition at line 208 of file GlobalLinSys.h.

References v_Initialise().

Referenced by Nektar::MultiRegions::GlobalLinSysStaticCond::v_InitObject(), and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_InitObject().

210  {
211  v_Initialise(pLocToGloMap);
212  }
virtual void v_Initialise(const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
void Nektar::MultiRegions::GlobalLinSys::InitObject ( )
inline

Definition at line 203 of file GlobalLinSys.h.

References v_InitObject().

204  {
205  v_InitObject();
206  }
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 180 of file GlobalLinSys.h.

References v_Solve().

185  {
186  v_Solve(in,out,locToGloMap,dirForcing);
187  }
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 193 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().

199  {
200  v_SolveLinearSystem(pNumRows, pInput, pOutput, locToGloMap, pNumDir);
201  }
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 368 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().

369  {
370  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
371 
372  StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(n);
373 
374  // need to be initialised with zero size for non variable
375  // coefficient case
376  StdRegions::VarCoeffMap vVarCoeffMap;
377 
378  // retrieve variable coefficients
379  if(m_linSysKey.GetNVarCoeffs() > 0)
380  {
381  StdRegions::VarCoeffMap::const_iterator x;
382  int cnt = expList->GetPhys_Offset(n);
383  for (x = m_linSysKey.GetVarCoeffs().begin();
384  x != m_linSysKey.GetVarCoeffs().end (); ++x)
385  {
386  vVarCoeffMap[x->first] = x->second + cnt;
387  }
388  }
389 
391  vExp->DetShapeType(),
392  *vExp,
394  vVarCoeffMap);
395 
396  vExp->DropLocStaticCondMatrix(matkey);
397  }
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:225
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:123
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
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 228 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().

229  {
230  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
231  int cnt = 0;
232  DNekScalMatSharedPtr loc_mat;
233 
235  boost::dynamic_pointer_cast<LocalRegions::Expansion>(
236  expList->GetExp(n));
237 
238  // need to be initialised with zero size for non variable
239  // coefficient case
240  StdRegions::VarCoeffMap vVarCoeffMap;
241 
242  // retrieve variable coefficients
243  if(m_linSysKey.GetNVarCoeffs() > 0)
244  {
245  StdRegions::VarCoeffMap::const_iterator x;
246  cnt = expList->GetPhys_Offset(n);
247 
248  for (x = m_linSysKey.GetVarCoeffs().begin();
249  x != m_linSysKey.GetVarCoeffs().end(); ++x)
250  {
251  vVarCoeffMap[x->first] = x->second + cnt;
252  }
253  }
254 
256  vExp->DetShapeType(),
257  *vExp, m_linSysKey.GetConstFactors(),
258  vVarCoeffMap);
259  loc_mat = vExp->GetLocMatrix(matkey);
260 
261  // apply robin boundary conditions to the matrix.
262  if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
263  {
265 
266  // declare local matrix from scaled matrix.
267  int rows = loc_mat->GetRows();
268  int cols = loc_mat->GetColumns();
269  const NekDouble *dat = loc_mat->GetRawPtr();
271  AllocateSharedPtr(rows,cols,dat);
272  Blas::Dscal(rows*cols,loc_mat->Scale(),new_mat->GetRawPtr(),1);
273 
274  // add local matrix contribution
275  for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
276  {
277  vExp->AddRobinMassMatrix(
278  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
279  }
280 
281  // redeclare loc_mat to point to new_mat plus the scalar.
283  1.0, new_mat);
284  }
285 
286  // finally return the matrix.
287  return loc_mat;
288  }
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:225
const map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:127
double NekDouble
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:123
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
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 216 of file GlobalLinSys.cpp.

References m_expList.

Referenced by GetNumBlocks().

217  {
218  return m_expList.lock()->GetExpSize();
219  }
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
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.

Definition at line 298 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().

300  {
301  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
302  int cnt = 0;
303  DNekScalBlkMatSharedPtr loc_mat;
304  DNekScalMatSharedPtr tmp_mat;
305 
306  StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(n);
307 
308  // need to be initialised with zero size for non variable
309  // coefficient case
310  StdRegions::VarCoeffMap vVarCoeffMap;
311 
312  // retrieve variable coefficients
313  if(m_linSysKey.GetNVarCoeffs() > 0)
314  {
315  StdRegions::VarCoeffMap::const_iterator x;
316  cnt = expList->GetPhys_Offset(n);
317  for (x = m_linSysKey.GetVarCoeffs().begin();
318  x != m_linSysKey.GetVarCoeffs().end (); ++x)
319  {
320  vVarCoeffMap[x->first] = x->second + cnt;
321  }
322  }
323 
325  vExp->DetShapeType(),
326  *vExp,
328  vVarCoeffMap);
329 
330  loc_mat = vExp->GetLocStaticCondMatrix(matkey);
331 
332  if(m_robinBCInfo.count(n) != 0) // add robin mass matrix
333  {
335 
336  tmp_mat = loc_mat->GetBlock(0,0);
337 
338  // declare local matrix from scaled matrix.
339  int rows = tmp_mat->GetRows();
340  int cols = tmp_mat->GetColumns();
341  const NekDouble *dat = tmp_mat->GetRawPtr();
343  AllocateSharedPtr(rows, cols, dat);
344  Blas::Dscal(rows*cols,tmp_mat->Scale(),new_mat->GetRawPtr(),1);
345 
346  // add local matrix contribution
347  for(rBC = m_robinBCInfo.find(n)->second;rBC; rBC = rBC->next)
348  {
349  vExp->AddRobinMassMatrix(
350  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, new_mat);
351  }
352 
353  // redeclare loc_mat to point to new_mat plus the scalar.
355  1.0, new_mat);
356  loc_mat->SetBlock(0,0,tmp_mat);
357  }
358 
359  return loc_mat;
360  }
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:225
const map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:127
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:123
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:125
void Nektar::MultiRegions::GlobalLinSys::v_Initialise ( const boost::shared_ptr< AssemblyMap > &  pLocToGloMap)
privatevirtual

Reimplemented in Nektar::MultiRegions::GlobalLinSysStaticCond.

Definition at line 404 of file GlobalLinSys.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by Initialise().

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

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

Definition at line 399 of file GlobalLinSys.cpp.

References ErrorUtil::efatal, and NEKERROR.

Referenced by InitObject().

400  {
401  NEKERROR(ErrorUtil::efatal, "Method does not exist" );
402  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:158
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 156 of file GlobalLinSys.h.

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

Definition at line 155 of file GlobalLinSys.h.

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

Local Matrix System.

Definition at line 125 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::GlobalLinSysIterativeFull::v_DoMatrixMultiply(), Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_DoMatrixMultiply(), v_DropStaticCondBlock(), v_GetBlock(), v_GetNumBlocks(), v_GetStaticCondBlock(), 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 map<int, RobinBCInfoSharedPtr> Nektar::MultiRegions::GlobalLinSys::m_robinBCInfo
protected

Robin boundary info.

Definition at line 127 of file GlobalLinSys.h.

Referenced by v_GetBlock(), and v_GetStaticCondBlock().