43 namespace MultiRegions
56 string GlobalLinSysIterativeFull::className
59 GlobalLinSysIterativeFull::create,
60 "Iterative solver for full matrix system.");
70 GlobalLinSysIterativeFull::GlobalLinSysIterativeFull(
72 const std::weak_ptr<ExpList> &pExp,
73 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
78 "This routine should only be used when using an Iterative " 79 "conjugate gradient matrix solve.");
116 std::shared_ptr<MultiRegions::ExpList> expList =
m_expList.lock();
120 if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
124 else if (std::dynamic_pointer_cast<AssemblyMapDG>(pLocToGloMap))
133 bool dirForcCalculated = (bool) pDirForcing.num_elements();
134 int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
135 int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
136 int nDirTotal = nDirDofs;
138 expList->GetComm()->GetRowComm()
146 if(dirForcCalculated)
149 pDirForcing.get(), 1,
156 expList->GeneralMatrixOp(
169 nGlobDofs, tmp, out, pLocToGloMap, nDirDofs);
170 Vmath::Vadd(nGlobDofs-nDirDofs, &out [nDirDofs], 1,
171 &pOutput[nDirDofs], 1, &pOutput[nDirDofs], 1);
175 ASSERTL0(
false,
"Need DG solve if using Dir BCs");
193 std::shared_ptr<MultiRegions::ExpList> expList =
m_expList.lock();
204 "Robin boundaries not set up in IterativeFull solver.");
205 int nGlobal = asmMap->GetNumGlobalCoeffs();
206 int nLocal = asmMap->GetNumLocalCoeffs();
207 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
208 int nNonDir = nGlobal - nDir;
217 asmMap->GlobalToLocal(pInput, robin_l);
221 for (
int n = 0; n < expList->GetNumElmts(); ++n)
224 int offset = expList->GetCoeff_Offset(n);
225 int ncoeffs = expList->GetExp(nel)->GetNcoeffs();
234 for(rBC =
m_robinBCInfo.find(nel)->second;rBC; rBC = rBC->next)
236 vExp->AddRobinEdgeContribution(rBC->m_robinID,rBC->m_robinPrimitiveCoeffs, tmp = robin_l + offset);
246 asmMap->LocalToGlobal(robin_l, robin_A);
248 Vmath::Vadd(nGlobal, pOutput, 1, robin_A, 1, pOutput, 1);
#define ASSERTL0(condition, msg)
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
virtual void v_UniqueMap()
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.
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
virtual ~GlobalLinSysIterativeFull()
Array< OneD, int > m_map
Global to universal unique map.
Describe a linear system.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
std::shared_ptr< Expansion > ExpansionSharedPtr
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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...
void Zero(int n, T *x, const int incx)
Zero vector.
GlobalLinSysFactory & GetGlobalLinSysFactory()
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
std::weak_ptr< AssemblyMap > m_locToGloMap
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)