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

#include <PreconditionerLinear.h>

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

Public Member Functions

 PreconditionerLinear (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
virtual ~PreconditionerLinear ()
 
- Public Member Functions inherited from Nektar::MultiRegions::Preconditioner
 Preconditioner (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
virtual ~Preconditioner ()
 
void DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce=NullNekDouble1DArray)
 
void DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
 
void DoTransformToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInOut)
 
void DoMultiplybyInverseTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoMultiplybyInverseTransposedTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void BuildPreconditioner ()
 
void InitObject ()
 
Array< OneD, NekDoubleAssembleStaticCondGlobalDiagonals ()
 Performs global assembly of diagonal entries to global Schur complement matrix. More...
 
const DNekScalBlkMatSharedPtrGetBlockTransformedSchurCompl () const
 
const DNekScalBlkMatSharedPtrGetBlockCMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockInvDMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockSchurCompl () const
 
const DNekScalBlkMatSharedPtrGetBlockTransformationMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockTransposedTransformationMatrix () const
 
DNekScalMatSharedPtr TransformedSchurCompl (int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
 

Static Public Member Functions

static PreconditionerSharedPtr create (const boost::shared_ptr< GlobalLinSys > &plinsys, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className1
 Name of class. More...
 

Protected Attributes

GlobalLinSysSharedPtr m_vertLinsys
 
boost::shared_ptr< AssemblyMapm_vertLocToGloMap
 
- Protected Attributes inherited from Nektar::MultiRegions::Preconditioner
const boost::weak_ptr
< GlobalLinSys
m_linsys
 
PreconditionerType m_preconType
 
DNekMatSharedPtr m_preconditioner
 
boost::shared_ptr< AssemblyMapm_locToGloMap
 
LibUtilities::CommSharedPtr m_comm
 

Private Member Functions

virtual void v_InitObject ()
 
virtual void v_DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
 Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom. More...
 
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Apply a preconditioner to the conjugate gradient method. More...
 
virtual void v_BuildPreconditioner ()
 

Static Private Attributes

static std::string solveType
 
static std::string solveTypeIds []
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::MultiRegions::Preconditioner
virtual DNekScalMatSharedPtr v_TransformedSchurCompl (int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
 Get block elemental transposed transformation matrix $\mathbf{R}^{T}$. More...
 

Detailed Description

This class implements preconditioning for the conjugate gradient matrix solver.

Definition at line 58 of file PreconditionerLinear.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::PreconditionerLinear::PreconditionerLinear ( const boost::shared_ptr< GlobalLinSys > &  plinsys,
const AssemblyMapSharedPtr pLocToGloMap 
)

Definition at line 86 of file PreconditionerLinear.cpp.

89  : Preconditioner(plinsys, pLocToGloMap)
90  {
91  }
Preconditioner(const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual Nektar::MultiRegions::PreconditionerLinear::~PreconditionerLinear ( )
inlinevirtual

Definition at line 80 of file PreconditionerLinear.h.

80 {}

Member Function Documentation

static PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerLinear::create ( const boost::shared_ptr< GlobalLinSys > &  plinsys,
const boost::shared_ptr< AssemblyMap > &  pLocToGloMap 
)
inlinestatic

Creates an instance of this class.

Definition at line 62 of file PreconditionerLinear.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

66  {
68  p->InitObject();
69  return p;
70  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:62
void Nektar::MultiRegions::PreconditionerLinear::v_BuildPreconditioner ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 97 of file PreconditionerLinear.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::eLinearPreconPETSc, Nektar::MultiRegions::eLinearPreconXxt, Nektar::StdRegions::eMass, Nektar::MultiRegions::ePETScFullMatrix, Nektar::MultiRegions::ePETScStaticCond, Nektar::StdRegions::ePreconLinearSpace, Nektar::StdRegions::ePreconLinearSpaceMass, Nektar::MultiRegions::eXxtFullMatrix, Nektar::MultiRegions::Preconditioner::m_linsys, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_vertLinsys, m_vertLocToGloMap, and solveType.

98  {
99  GlobalSysSolnType sType = m_locToGloMap->GetGlobalSysSolnType();
100  ASSERTL0(sType == eIterativeStaticCond || sType == ePETScStaticCond,
101  "This type of preconditioning is not implemented "
102  "for this solver");
103 
104  boost::shared_ptr<MultiRegions::ExpList>
105  expList=((m_linsys.lock())->GetLocMat()).lock();
106 
108  expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
109  "LinearPreconSolver");
110 
111  GlobalSysSolnType linSolveType;
112 
113  switch(solveType)
114  {
115  case eLinearPreconPETSc:
116  {
117  linSolveType = ePETScFullMatrix;
118 #ifndef NEKTAR_USING_PETSC
119  ASSERTL0(false, "Nektar++ has not been compiled with "
120  "PETSc support.");
121 #endif
122  }
123  case eLinearPreconXxt:
124  default:
125  {
126  linSolveType = eXxtFullMatrix;
127  break;
128  }
129  }
130 
131  m_vertLocToGloMap = m_locToGloMap->LinearSpaceMap(*expList, linSolveType);
132 
133  // Generate linear solve system.
134  StdRegions::MatrixType mType =
135  m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass ?
138 
139  GlobalLinSysKey preconKey(mType, m_vertLocToGloMap, (m_linsys.lock())->GetKey().GetConstFactors());
140 
141  switch(solveType)
142  {
143  case eLinearPreconXxt:
144  {
146  AllocateSharedPtr(preconKey,expList,m_vertLocToGloMap);
147  break;
148  }
149  case eLinearPreconPETSc:
150  {
151 #ifdef NEKTAR_USING_PETSC
153  AllocateSharedPtr(preconKey,expList,m_vertLocToGloMap);
154 #else
155  ASSERTL0(false, "Nektar++ has not been compiled with "
156  "PETSc support.");
157 #endif
158  }
159  }
160  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< AssemblyMap > m_vertLocToGloMap
const boost::weak_ptr< GlobalLinSys > m_linsys
void Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Apply a preconditioner to the conjugate gradient method.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 165 of file PreconditionerLinear.cpp.

References Nektar::NullNekDouble1DArray, and v_DoPreconditionerWithNonVertOutput().

168  {
171  }
static Array< OneD, NekDouble > NullNekDouble1DArray
virtual void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...
void Nektar::MultiRegions::PreconditionerLinear::v_DoPreconditionerWithNonVertOutput ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const Array< OneD, NekDouble > &  pNonVertOutput,
Array< OneD, NekDouble > &  pVertForce 
)
privatevirtual

Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of freedom.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 176 of file PreconditionerLinear.cpp.

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScStaticCond, Nektar::MultiRegions::Preconditioner::m_locToGloMap, m_vertLinsys, m_vertLocToGloMap, Nektar::NullNekDouble1DArray, Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_DoPreconditioner().

181  {
182  GlobalSysSolnType solvertype=m_locToGloMap->GetGlobalSysSolnType();
183  switch(solvertype)
184  {
187  {
188  int i,val;
189  int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
190  int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
191  // mapping from full space to vertices
192  Array<OneD, int> LocToGloBnd = m_vertLocToGloMap->GetLocalToGlobalBndMap();
193 
194  // Global to local for linear solver (different from above)
195  Array<OneD, int> LocToGlo = m_vertLocToGloMap->GetLocalToGlobalMap();
196 
197  // number of Dir coeffs in from full problem
198  int nDirFull = m_locToGloMap->GetNumGlobalDirBndCoeffs();
199 
200  Array<OneD,NekDouble> In(nglo,0.0);
201  Array<OneD,NekDouble> Out(nglo,0.0);
202 
203  // Gather rhs
204  for(i = 0; i < nloc; ++i)
205  {
206  val = LocToGloBnd[i];
207  if(val >= nDirFull)
208  {
209  In[LocToGlo[i]] = pInput[val-nDirFull];
210  }
211  }
212 
213  // Do solve without enforcing any boundary conditions.
214  m_vertLinsys->SolveLinearSystem(
215  m_vertLocToGloMap->GetNumGlobalCoeffs(),
216  In,Out,m_vertLocToGloMap,
217  m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
218 
219 
220  if(pNonVertOutput != NullNekDouble1DArray)
221  {
222  ASSERTL1(pNonVertOutput.num_elements() >= pOutput.num_elements(),"Non Vert output is not of sufficient length");
223  Vmath::Vcopy(pOutput.num_elements(),pNonVertOutput,1,pOutput,1);
224  }
225  else
226  {
227  //Copy input to output as a unit preconditioner on
228  //any other value
229  Vmath::Vcopy(pInput.num_elements(),pInput,1,pOutput,1);
230  }
231 
232  if(pVertForce != NullNekDouble1DArray)
233  {
234  Vmath::Zero(pVertForce.num_elements(),pVertForce,1);
235  // Scatter back soln from linear solve
236  for(i = 0; i < nloc; ++i)
237  {
238  val = LocToGloBnd[i];
239  if(val >= nDirFull)
240  {
241  pOutput[val-nDirFull] = Out[LocToGlo[i]];
242  // copy vertex forcing into this vector
243  pVertForce[val-nDirFull] = In[LocToGlo[i]];
244  }
245  }
246  }
247  else
248  {
249  // Scatter back soln from linear solve
250  for(i = 0; i < nloc; ++i)
251  {
252  val = LocToGloBnd[i];
253  if(val >= nDirFull)
254  {
255  pOutput[val-nDirFull] = Out[LocToGlo[i]];
256  }
257  }
258  }
259  }
260  break;
261  default:
262  ASSERTL0(0,"Unsupported solver type");
263  break;
264  }
265  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
static Array< OneD, NekDouble > NullNekDouble1DArray
boost::shared_ptr< AssemblyMap > m_vertLocToGloMap
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::MultiRegions::PreconditionerLinear::v_InitObject ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 93 of file PreconditionerLinear.cpp.

94  {
95  }

Member Data Documentation

string Nektar::MultiRegions::PreconditionerLinear::className1
static
Initial value:
"FullLinearSpace",
"Full Linear space inverse Preconditioning")

Name of class.

Registers the class with the Factory.

Definition at line 73 of file PreconditionerLinear.h.

GlobalLinSysSharedPtr Nektar::MultiRegions::PreconditionerLinear::m_vertLinsys
protected
boost::shared_ptr<AssemblyMap> Nektar::MultiRegions::PreconditionerLinear::m_vertLocToGloMap
protected
std::string Nektar::MultiRegions::PreconditionerLinear::solveType
staticprivate
Initial value:

Definition at line 87 of file PreconditionerLinear.h.

Referenced by v_BuildPreconditioner().

std::string Nektar::MultiRegions::PreconditionerLinear::solveTypeIds
staticprivate
Initial value:

Definition at line 88 of file PreconditionerLinear.h.