38 #include <boost/core/ignore_unused.hpp>    49 string DriverModifiedArnoldi::className =
    51                                     DriverModifiedArnoldi::create);
    52 string DriverModifiedArnoldi::driverLookupId =
    53         LibUtilities::SessionReader::RegisterEnumValue(
"Driver",
    59 DriverModifiedArnoldi::DriverModifiedArnoldi(
    82     m_equ[0]->PrintSummary(out);
    85     if (
m_comm->GetRank() == 0)
    87         out << 
"\tArnoldi solver type    : Modified Arnoldi" << endl;
    92     for( 
int i = 0; i < 
m_nequ; ++i)
    94         m_equ[i]->DoInitialise();
    98     m_equ[m_nequ-1] ->TransPhysToCoeff();
   110     int nq              = 
m_equ[0]->UpdateFields()[0]->GetNcoeffs();
   115     std::string evlFile = 
m_session->GetSessionName() + 
".evl";
   117     if (
m_comm->GetRank() == 0)
   119         evlout.open(evlFile.c_str());
   132     for (i = 0; i < m_kdim + 1; ++i)
   139     if(
m_session->DefinesFunction(
"InitialConditions"))
   141         if (
m_comm->GetRank() == 0)
   143             out << 
"\tInital vector       : specified in input file " << endl;
   145         m_equ[0]->SetInitialConditions(0.0,
false);
   151         if (
m_comm->GetRank() == 0)
   153             out << 
"\tInital vector       : random  " << endl;
   163     if (
m_comm->GetRank() == 0)
   165         out << 
"Iteration: " << 0 <<  endl;
   169     alpha[0] = 
Blas::Ddot(ntot, &Kseq[0][0], 1, &Kseq[0][0], 1);
   171     alpha[0] = std::sqrt(alpha[0]);
   172     Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
   176     for (i = 1; !converged && i <= 
m_kdim; ++i)
   182         alpha[i] = 
Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
   184         alpha[i] = std::sqrt(alpha[i]);
   187         Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
   190         for (
int k = 0; k < i + 1; ++k)
   196         EV_small(Tseq, ntot, alpha, i, zvec, wr, wi, resnorm);
   199         converged = 
EV_test(i, i, zvec, wr, wi, resnorm,
   200             std::min(i, 
m_nvec), evlout, resid0);
   203             converged = max (converged, 0);
   210         if (
m_comm->GetRank() == 0)
   212             out << 
"Iteration: " <<  i << 
" (residual : " << resid0
   220         for (i = m_kdim + 1; !converged && i <= 
m_nits; ++i)
   224             for (
int j = 1; j <= 
m_kdim; ++j)
   226                 alpha[j-1] = alpha[j];
   231             EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
   235                                              &Kseq[m_kdim][0], 1);
   237             alpha[
m_kdim] = std::sqrt(alpha[m_kdim]);
   238             Vmath::Smul(ntot, 1.0/alpha[m_kdim], Kseq[m_kdim], 1,
   242             for (
int k = 0; k < m_kdim + 1; ++k)
   248             EV_small(Tseq, ntot, alpha, m_kdim, zvec, wr, wi, resnorm);
   251             converged = 
EV_test(i, m_kdim, zvec, wr, wi, resnorm,
   254             if (
m_comm->GetRank() == 0)
   256                 out << 
"Iteration: " <<  i << 
" (residual : "   257                     << resid0 << 
")" <<endl;
   268     for(j = 0; j < 
m_equ[0]->GetNvariables(); ++j)
   272         if (
m_comm->GetRank() == 0)
   274             out << 
"L 2 error (variable " << 
m_equ[0]->GetVariable(j)
   275             << 
") : " << vL2Error << endl;
   276             out << 
"L inf error (variable " << 
m_equ[0]->GetVariable(j)
   277             << 
") : " << vLinfError << endl;
   282     EV_post(Tseq, Kseq, ntot, min(--i, m_kdim), 
m_nvec, zvec, wr, wi,
   286               "Need to implement Ritz re-evaluation of"   287               "eigenvalue. Only one half of complex "   288               "value will be correct");
   295     if (
m_comm->GetRank() == 0)
   311     m_equ[0]->TransCoeffToPhys();
   313     m_equ[0]->SetTime(0.);
   319         fields = 
m_equ[0]->UpdateFields();
   323         m_equ[1]->TransCoeffToPhys();
   325         m_equ[1]->SetTime(0.);
   348     int kdimp = kdim + 1;
   356     for (
int i = 0; i < kdimp; ++i)
   360         gsc = std::sqrt(gsc);
   361         ASSERTL0(gsc != 0.0, 
"Vectors are linearly independent.");
   364         Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
   366         for (
int j = i + 1; j < kdimp; ++j)
   368             gsc = 
Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
   370             Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
   376     for (
int i = 0; i < kdim; ++i)
   378         for (
int j = 0; j < kdim; ++j)
   380             H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
   381             - 
Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
   382             H[j*kdim+i] /= R[j*kdimp+j];
   386     H[(kdim-1)*kdim+kdim] = alpha[kdim]
   387                   * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
   389     Lapack::dgeev_(
'N', 
'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
   390                    &zvec[0], kdim, &rwork[0], lwork, ier);
   394     resnorm = H[(kdim-1)*kdim + kdim];
   416     for (
int i = 0; i < kdim; ++i)
   419                                              &zvec[0] + i*kdim, 1));
   420         resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
   423             resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
   427     EV_sort(zvec, wr, wi, resid, kdim);
   434     if (
m_comm->GetRank() == 0)
   436         evlout << 
"-- Iteration = " << itrn << 
", H(k+1, k) = "   439         evlout.setf(ios::scientific, ios::floatfield);
   442             evlout << 
"        Magnitude   Angle       Growth      "   443                    << 
"Frequency   Residual" << endl;
   447             evlout << 
"        Real        Imaginary   inverse real  "   448                    << 
"inverse imag  Residual" << endl;
   451         for (
int i = 0; i < kdim; i++)
   453             WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
   457     resid0 = resid[nvec-1];
   474     for (
int j = 1; j < dim; ++j)
   481         while (i >= 0 && test[i] > te_tmp)
   486             Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
   492         Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
   516         ASSERTL0(
false, 
"Convergence was not achieved within the "   517                         "prescribed number of iterations.");
   522         ASSERTL0(
false, 
"Minimum residual reached.");
   524     else if (icon == nvec)
   527         EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
   529                 = 
m_equ[0]->UpdateFields();
   531         for (
int j = 0; j < icon; ++j)
   533             std::string file = 
m_session->GetSessionName() + 
"_eig_"   534                 + boost::lexical_cast<std::string>(j)
   537             if (
m_comm->GetRank() == 0)
   547         ASSERTL0(
false, 
"Unrecognised value.");
   565     boost::ignore_unused(wr);
   570     for (
int j = 0; j < nvec; ++j)
   573         for (
int i = 0; i < kdim; ++i)
   575             wgt = zvec[i + j*kdim];
   576             Vmath::Svtvp(ntot, wgt, bvecs[i], 1, evecs[j], 1, evecs[j], 1);
   581     for (
int i = 0; i < nvec; ++i)
   585             norm = 
Blas::Ddot(ntot, &evecs[i][0], 1, &evecs[i][0], 1);
   587             norm = std::sqrt(norm);
   588             Vmath::Smul(ntot, 1.0/norm, evecs[i], 1, evecs[i], 1);
   592             norm  = 
Blas::Ddot(ntot, &evecs[i][0],   1, &evecs[i][0],   1);
   593             norm += 
Blas::Ddot(ntot, &evecs[i+1][0], 1, &evecs[i+1][0], 1);
   595             norm = std::sqrt(norm);
   597             Vmath::Smul(ntot, 1.0/norm, evecs[i],   1, evecs[i],   1);
   598             Vmath::Smul(ntot, 1.0/norm, evecs[i+1], 1, evecs[i+1], 1);
 Array< OneD, NekDouble > m_imag_evl
#define ASSERTL0(condition, msg)
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations. 
NekDouble m_evtol
Maxmum number of iterations. 
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields. 
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, std::ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H. 
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y 
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
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. 
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage. 
virtual void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation. 
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file. 
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. 
virtual ~DriverModifiedArnoldi()
Destructor. 
void EV_update(Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
Generates a new vector in the sequence by applying the linear operator. 
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y. 
int m_nvec
Dimension of Krylov subspace. 
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated. 
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
#define WARNINGL0(condition, msg)
virtual void v_InitObject(std::ostream &out=std::cout)
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator. 
int m_nfields
interval to dump information if required. 
bool m_timeSteppingAlgorithm
Period of time stepping algorithm. 
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y 
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = . 
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve. 
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)
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise. 
DriverFactory & GetDriverFactory()
LibUtilities::CommSharedPtr m_comm
Communication object. 
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory. 
virtual void v_Execute(std::ostream &out=std::cout)
Virtual function for solve implementation. 
Base class for the development of solvers. 
void Zero(int n, T *x, const int incx)
Zero vector. 
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
LibUtilities::SessionReaderSharedPtr m_session
Session reader object. 
int m_nequ
number of equations 
int m_nits
Number of vectors to test.