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

#include <DriverModifiedArnoldi.h>

Inheritance diagram for Nektar::SolverUtils::DriverModifiedArnoldi:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::DriverModifiedArnoldi:
Collaboration graph
[legend]

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 Creates an instance of this class.

Static Public Attributes

static std::string className = GetDriverFactory().RegisterCreatorFunction("ModifiedArnoldi", DriverModifiedArnoldi::create)
 Name of the class.

Protected Member Functions

 DriverModifiedArnoldi (const LibUtilities::SessionReaderSharedPtr pSession)
 Constructor.
virtual ~DriverModifiedArnoldi ()
 Destructor.
virtual void v_InitObject (ostream &out=cout)
 Virtual function for initialisation implementation.
virtual void v_Execute (ostream &out=cout)
 Virtual function for solve implementation.
- Protected Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
 DriverArnoldi (const LibUtilities::SessionReaderSharedPtr pSession)
 Constructor.
virtual ~DriverArnoldi ()
 Destructor.
void CopyArnoldiArrayToField (Array< OneD, NekDouble > &array)
 Copy Arnoldi storage to fields.
void CopyFieldToArnoldiArray (Array< OneD, NekDouble > &array)
 Copy fields to Arnoldi storage.
void CopyFwdToAdj ()
 Copy the forward field to the adjoint system in transient growth calculations.
void WriteFld (std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
void WriteFld (std::string file, Array< OneD, NekDouble > coeffs)
void WriteEvs (ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble)
virtual Array< OneD, NekDoublev_GetRealEvl (void)
virtual Array< OneD, NekDoublev_GetImagEvl (void)
- Protected Member Functions inherited from Nektar::SolverUtils::Driver
 Driver (const LibUtilities::SessionReaderSharedPtr pSession)
 Initialises EquationSystem class members.

Private Member Functions

void EV_update (Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
 Generates a new vector in the sequence by applying the linear operator.
void EV_small (Array< OneD, Array< OneD, NekDouble > > &Kseq, const int ntot, const Array< OneD, NekDouble > &alpha, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, NekDouble &resnorm)
 Generates the upper Hessenberg matrix H and computes its eigenvalues.
int EV_test (const int itrn, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const NekDouble resnorm, const int nvec, ofstream &evlout, NekDouble &resid0)
 Tests for convergence of eigenvalues of H.
void EV_sort (Array< OneD, NekDouble > &evec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, Array< OneD, NekDouble > &test, const int dim)
 Sorts a sequence of eigenvectors/eigenvalues by magnitude.
void EV_post (Array< OneD, Array< OneD, NekDouble > > &Tseq, Array< OneD, Array< OneD, NekDouble > > &Kseq, const int ntot, const int kdim, const int nvec, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const int icon)
void EV_big (Array< OneD, Array< OneD, NekDouble > > &bvecs, Array< OneD, Array< OneD, NekDouble > > &evecs, const int ntot, const int kdim, const int nvec, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi)

Static Private Attributes

static std::string driverLookupId = LibUtilities::SessionReader::RegisterEnumValue("Driver","ModifiedArnoldi",0)

Friends

class MemoryManager< DriverModifiedArnoldi >

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
SOLVER_UTILS_EXPORT void ArnoldiSummary (std::ostream &out)
- Protected Attributes inherited from Nektar::SolverUtils::DriverArnoldi
int m_kdim
int m_nvec
 Dimension of Krylov subspace.
int m_nits
 Number of vectors to test.
NekDouble m_evtol
 Maxmum number of iterations.
NekDouble m_period
 Tolerance of iteratiosn.
bool m_timeSteppingAlgorithm
 Period of time stepping algorithm.
int m_infosteps
 underlying operator is time stepping
int m_nfields
 interval to dump information if required.
NekDouble m_realShift
NekDouble m_imagShift
Array< OneD, NekDoublem_real_evl
Array< OneD, NekDoublem_imag_evl
- Static Protected Attributes inherited from Nektar::SolverUtils::Driver
static std::string evolutionOperatorLookupIds []
static std::string evolutionOperatorDef = LibUtilities::SessionReader::RegisterDefaultSolverInfo("EvolutionOperator","Nonlinear")
static std::string driverDefault = LibUtilities::SessionReader::RegisterDefaultSolverInfo("Driver","Standard")

Detailed Description

Definition at line 46 of file DriverModifiedArnoldi.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::DriverModifiedArnoldi::DriverModifiedArnoldi ( const LibUtilities::SessionReaderSharedPtr  pSession)
protected

Constructor.

Definition at line 51 of file DriverModifiedArnoldi.cpp.

: DriverArnoldi(pSession)
{
}
Nektar::SolverUtils::DriverModifiedArnoldi::~DriverModifiedArnoldi ( )
protectedvirtual

Destructor.

Definition at line 61 of file DriverModifiedArnoldi.cpp.

{
}

Member Function Documentation

static DriverSharedPtr Nektar::SolverUtils::DriverModifiedArnoldi::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file DriverModifiedArnoldi.h.

{
p->InitObject();
return p;
}
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_big ( Array< OneD, Array< OneD, NekDouble > > &  bvecs,
Array< OneD, Array< OneD, NekDouble > > &  evecs,
const int  ntot,
const int  kdim,
const int  nvec,
Array< OneD, NekDouble > &  zvec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi 
)
private

Definition at line 476 of file DriverModifiedArnoldi.cpp.

References Vmath::Dot(), Nektar::Dot(), Vmath::Smul(), Vmath::Svtvp(), and Vmath::Zero().

Referenced by EV_post().

{
NekDouble wgt, norm;
// Generate big eigenvectors
for (int j = 0; j < nvec; ++j)
{
Vmath::Zero(ntot, evecs[j], 1);
for (int i = 0; i < kdim; ++i)
{
wgt = zvec[i + j*kdim];
Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
}
}
// Normalise the big eigenvectors
for (int i = 0; i < nvec; ++i)
{
if (wi[i] == 0.0) // Real mode
{
norm = std::sqrt(Vmath::Dot(ntot, evecs[i], evecs[i]));
Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
}
else
{
norm = Vmath::Dot(ntot, evecs[i], 1, evecs[i], 1);
norm += Vmath::Dot(ntot, evecs[i+1], 1, evecs[i+1], 1);
norm = std::sqrt(norm);
Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
i++;
}
}
}
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_post ( Array< OneD, Array< OneD, NekDouble > > &  Tseq,
Array< OneD, Array< OneD, NekDouble > > &  Kseq,
const int  ntot,
const int  kdim,
const int  nvec,
Array< OneD, NekDouble > &  zvec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi,
const int  icon 
)
private

Definition at line 430 of file DriverModifiedArnoldi.cpp.

References ASSERTL0, EV_big(), Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::Driver::m_session, and Nektar::SolverUtils::DriverArnoldi::WriteFld().

Referenced by v_Execute().

{
if (icon == 0)
{
// Not converged, write final Krylov vector
ASSERTL0(false, "Convergence was not achieved within the prescribed number of iterations.");
}
else if (icon < 0)
{
// Minimum residual reached
ASSERTL0(false, "Minimum residual reached.");
}
else if (icon == nvec)
{
// Converged, write out eigenvectors
EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
Array<OneD, MultiRegions::ExpListSharedPtr> fields
= m_equ[0]->UpdateFields();
for (int j = 0; j < icon; ++j)
{
std::string file = m_session->GetFilename().substr(0,m_session->GetFilename().find_last_of('.')) + "_eig_" + boost::lexical_cast<std::string>(j);
WriteFld(file,Kseq[j]);
}
}
else
{
// Not recognised
ASSERTL0(false, "Unrecognised value.");
}
}
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_small ( Array< OneD, Array< OneD, NekDouble > > &  Kseq,
const int  ntot,
const Array< OneD, NekDouble > &  alpha,
const int  kdim,
Array< OneD, NekDouble > &  zvec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi,
NekDouble resnorm 
)
private

Generates the upper Hessenberg matrix H and computes its eigenvalues.

Definition at line 284 of file DriverModifiedArnoldi.cpp.

References ASSERTL0, Vmath::Dot(), Nektar::Dot(), Vmath::Smul(), and Vmath::Svtvp().

Referenced by v_Execute().

{
int kdimp = kdim + 1;
int lwork = 10*kdim;
int ier;
Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
Array<OneD, NekDouble> rwork(lwork, 0.0);
// Modified G-S orthonormalisation
for (int i = 0; i < kdimp; ++i)
{
NekDouble gsc = std::sqrt(Vmath::Dot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1));
ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
R[i*kdimp+i] = gsc;
Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
for (int j = i + 1; j < kdimp; ++j)
{
gsc = Vmath::Dot(ntot, Kseq[i], 1, Kseq[j], 1);
Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
R[j*kdimp+i] = gsc;
}
}
// Compute H matrix
for (int i = 0; i < kdim; ++i)
{
for (int j = 0; j < kdim; ++j)
{
H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
- Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
H[j*kdim+i] /= R[j*kdimp+j];
}
}
H[(kdim-1)*kdim+kdim] = alpha[kdim]
* std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
Lapack::dgeev_('N','V',kdim,&H[0],kdim,&wr[0],&wi[0],0,1,&zvec[0],kdim,&rwork[0],lwork,ier);
ASSERTL0(!ier, "Error with dgeev");
resnorm = H[(kdim-1)*kdim + kdim];
}
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_sort ( Array< OneD, NekDouble > &  evec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi,
Array< OneD, NekDouble > &  test,
const int  dim 
)
private

Sorts a sequence of eigenvectors/eigenvalues by magnitude.

Definition at line 395 of file DriverModifiedArnoldi.cpp.

References Vmath::Vcopy().

Referenced by EV_test().

{
Array<OneD, NekDouble> z_tmp(dim,0.0);
NekDouble wr_tmp, wi_tmp, te_tmp;
for (int j = 1; j < dim; ++j)
{
wr_tmp = wr[j];
wi_tmp = wi[j];
te_tmp = test[j];
Vmath::Vcopy(dim, &evec[0] + j*dim, 1, &z_tmp[0], 1);
int i = j - 1;
while (i >= 0 && test[i] > te_tmp)
{
wr[i+1] = wr[i];
wi[i+1] = wi[i];
test[i+1] = test[i];
Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
i--;
}
wr[i+1] = wr_tmp;
wi[i+1] = wi_tmp;
test[i+1] = te_tmp;
Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
}
};
int Nektar::SolverUtils::DriverModifiedArnoldi::EV_test ( const int  itrn,
const int  kdim,
Array< OneD, NekDouble > &  zvec,
Array< OneD, NekDouble > &  wr,
Array< OneD, NekDouble > &  wi,
const NekDouble  resnorm,
const int  nvec,
ofstream &  evlout,
NekDouble resid0 
)
private

Tests for convergence of eigenvalues of H.

Definition at line 341 of file DriverModifiedArnoldi.cpp.

References Vmath::Dot(), EV_sort(), Nektar::SolverUtils::DriverArnoldi::m_evtol, Nektar::SolverUtils::DriverArnoldi::m_timeSteppingAlgorithm, and Nektar::SolverUtils::DriverArnoldi::WriteEvs().

Referenced by v_Execute().

{
int idone = 0;
// NekDouble period = 0.1;
Array<OneD, NekDouble> resid(kdim);
for (int i = 0; i < kdim; ++i)
{
resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) /
std::sqrt(Vmath::Dot(kdim, &zvec[0] + i*kdim, 1, &zvec[0] + i*kdim, 1));
if (wi[i] < 0.0) resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
}
EV_sort(zvec, wr, wi, resid, kdim);
if (resid[nvec-1] < m_evtol) idone = nvec;
evlout << "-- Iteration = " << itrn << ", H(k+1, k) = " << resnorm << endl;
evlout.precision(4);
evlout.setf(ios::scientific, ios::floatfield);
{
evlout << "EV Magnitude Angle Growth Frequency Residual"
<< endl;
}
else
{
evlout << "EV Real Imaginary inverse real inverse imag Residual"
<< endl;
}
for (int i = 0; i < kdim; i++)
{
WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
}
resid0 = resid[nvec-1];
return idone;
}
void Nektar::SolverUtils::DriverModifiedArnoldi::EV_update ( Array< OneD, NekDouble > &  src,
Array< OneD, NekDouble > &  tgt 
)
private

Generates a new vector in the sequence by applying the linear operator.

Definition at line 251 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::CopyArnoldiArrayToField(), Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Nektar::SolverUtils::DriverArnoldi::CopyFwdToAdj(), Nektar::SolverUtils::eTransientGrowth, Nektar::SolverUtils::Driver::m_equ, and Nektar::SolverUtils::Driver::m_EvolutionOperator.

Referenced by v_Execute().

{
// Copy starting vector into first sequence element.
m_equ[0]->TransCoeffToPhys();
m_equ[0]->DoSolve();
{
Array<OneD, MultiRegions::ExpListSharedPtr> fields;
fields = m_equ[0]->UpdateFields();
//start Adjoint with latest fields of direct
m_equ[1]->TransCoeffToPhys();
m_equ[1]->DoSolve();
}
// Copy starting vector into first sequence element.
}
void Nektar::SolverUtils::DriverModifiedArnoldi::v_Execute ( ostream &  out = cout)
protectedvirtual

Virtual function for solve implementation.

Implements Nektar::SolverUtils::Driver.

Definition at line 91 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Vmath::Ddot(), Vmath::Dot(), EV_post(), EV_small(), EV_test(), EV_update(), Vmath::FillWhiteNoise(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_equ, Nektar::SolverUtils::DriverArnoldi::m_imag_evl, Nektar::SolverUtils::DriverArnoldi::m_kdim, Nektar::SolverUtils::DriverArnoldi::m_nfields, Nektar::SolverUtils::DriverArnoldi::m_nits, Nektar::SolverUtils::DriverArnoldi::m_nvec, Nektar::SolverUtils::DriverArnoldi::m_real_evl, Nektar::SolverUtils::Driver::m_session, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), and Vmath::Vcopy().

{
int i = 0;
int j = 0;
int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
int ntot = m_nfields*nq;
int converged = 0;
NekDouble resnorm = 0.0;
std::string evlFile = m_session->GetFilename().substr(0,m_session->GetFilename().find_last_of('.')) + ".evl";
ofstream evlout(evlFile.c_str());
// Allocate memory
Array<OneD, NekDouble> alpha = Array<OneD, NekDouble> (m_kdim+1, 0.0);
Array<OneD, NekDouble> wr = Array<OneD, NekDouble> (m_kdim, 0.0);
Array<OneD, NekDouble> wi = Array<OneD, NekDouble> (m_kdim, 0.0);
Array<OneD, NekDouble> zvec = Array<OneD, NekDouble> (m_kdim*m_kdim, 0.0);
Array<OneD, Array<OneD, NekDouble> > Kseq
= Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
Array<OneD, Array<OneD, NekDouble> > Tseq
= Array<OneD, Array<OneD, NekDouble> > (m_kdim + 1);
for (i = 0; i < m_kdim + 1; ++i)
{
Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
}
// Copy starting vector into second sequence element (temporary).
if(m_session->DefinesFunction("InitialConditions"))
{
out << "\tInital vector : specified in input file " << endl;
m_equ[0]->SetInitialConditions(0.0,false);
}
else
{
out << "\tInital vector : random " << endl;
NekDouble eps=1;
Vmath::FillWhiteNoise(ntot, eps , &Kseq[1][0], 1);
}
// Perform one iteration to enforce boundary conditions.
// Set this as the initial value in the sequence.
EV_update(Kseq[1], Kseq[0]);
out << "Iteration: " << 0 << endl;
// Normalise first vector in sequence
alpha[0] = std::sqrt(Blas::Ddot(ntot, &Kseq[0][0], 1,
&Kseq[0][0], 1));
m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
//alpha[0] = std::sqrt(alpha[0]);
Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
// Fill initial krylov sequence
NekDouble resid0;
for (i = 1; !converged && i <= m_kdim; ++i)
{
// Compute next vector
EV_update(Kseq[i-1], Kseq[i]);
// Normalise
alpha[i] = std::sqrt(Blas::Ddot(ntot, &Kseq[i][0], 1,
&Kseq[i][0], 1));
m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
//alpha[i] = std::sqrt(alpha[i]);
Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
// Copy Krylov sequence into temporary storage
for (int k = 0; k < i + 1; ++k)
{
Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
}
// Generate Hessenberg matrix and compute eigenvalues of it.
EV_small(Tseq, ntot, alpha, i, zvec, wr, wi, resnorm);
// Test for convergence.
converged = EV_test(i,i,zvec,wr,wi,resnorm,std::min(i,m_nvec),evlout,resid0);
converged = max (converged, 0);
out << "Iteration: " << i << " (residual : " << resid0 << ")" <<endl;
}
// Continue with full sequence
if (!converged)
{
for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
{
// Shift all the vectors in the sequence.
// First vector is removed.
//NekDouble invnorm = 1.0/sqrt(Blas::Ddot(ntot,Kseq[1],1,Kseq[1],1));
for (int j = 1; j <= m_kdim; ++j)
{
alpha[j-1] = alpha[j];
//Vmath::Smul(ntot,invnorm,Kseq[j],1,Kseq[j],1);
Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j-1], 1);
}
// Compute next vector
EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
// Compute new scale factor
alpha[m_kdim] = std::sqrt(Vmath::Dot(ntot, &Kseq[m_kdim][0], 1, &Kseq[m_kdim][0], 1));
//alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1, Kseq[m_kdim], 1);
// Copy Krylov sequence into temporary storage
for (int k = 0; k < m_kdim + 1; ++k)
{
Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
}
// Generate Hessenberg matrix and compute eigenvalues of it
EV_small(Tseq, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
// Test for convergence.
converged = EV_test(i,m_kdim,zvec,wr,wi,resnorm,m_nvec,evlout,resid0);
out << "Iteration: " << i << " (residual : " << resid0 << ")" <<endl;
}
}
m_equ[0]->Output();
// Evaluate and output computation time and solution accuracy.
// The specific format of the error output is essential for the
// regression tests to work.
// Evaluate L2 Error
for(j = 0; j < m_equ[0]->GetNvariables(); ++j)
{
NekDouble vL2Error = m_equ[0]->L2Error(j,false);
NekDouble vLinfError = m_equ[0]->LinfError(j);
if (m_comm->GetRank() == 0)
{
out << "L 2 error (variable " << m_equ[0]->GetVariable(j) << ") : " << vL2Error << endl;
out << "L inf error (variable " << m_equ[0]->GetVariable(j) << ") : " << vLinfError << endl;
}
}
// Process eigenvectors and write out.
EV_post(Tseq, Kseq, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi, converged);
// store eigenvalues so they can be access from driver class
m_real_evl = wr;
m_imag_evl = wi;
// Close the runtime info file.
evlout.close();
}
void Nektar::SolverUtils::DriverModifiedArnoldi::v_InitObject ( ostream &  out = cout)
protectedvirtual

Virtual function for initialisation implementation.

Reimplemented from Nektar::SolverUtils::DriverArnoldi.

Definition at line 69 of file DriverModifiedArnoldi.cpp.

References Nektar::SolverUtils::DriverArnoldi::ArnoldiSummary(), Nektar::SolverUtils::Driver::m_equ, and Nektar::SolverUtils::Driver::m_nequ.

{
m_equ[0]->PrintSummary(out);
// Print session parameters
out << "\tArnoldi solver type : Modified Arnoldi" << endl;
m_equ[m_nequ - 1]->DoInitialise();
//FwdTrans Initial conditions to be in Coefficient Space
m_equ[m_nequ-1] ->TransPhysToCoeff();
}

Friends And Related Function Documentation

friend class MemoryManager< DriverModifiedArnoldi >
friend

Definition at line 49 of file DriverModifiedArnoldi.h.

Member Data Documentation

string Nektar::SolverUtils::DriverModifiedArnoldi::className = GetDriverFactory().RegisterCreatorFunction("ModifiedArnoldi", DriverModifiedArnoldi::create)
static

Name of the class.

Definition at line 59 of file DriverModifiedArnoldi.h.

string Nektar::SolverUtils::DriverModifiedArnoldi::driverLookupId = LibUtilities::SessionReader::RegisterEnumValue("Driver","ModifiedArnoldi",0)
staticprivate

Definition at line 128 of file DriverModifiedArnoldi.h.