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();
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());
134 for (i = 0; i <
m_kdim + 1; ++i)
143 Kseqcopy[i] = Kseq[i];
149 if(
m_session->DefinesFunction(
"InitialConditions"))
151 if (
m_comm->GetRank() == 0)
153 out <<
"\tInital vector : specified in input file " << endl;
155 m_equ[0]->SetInitialConditions(0.0,
false);
161 if (
m_comm->GetRank() == 0)
163 out <<
"\tInital vector : random " << endl;
177 if (
m_comm->GetRank() == 0)
179 out <<
"Iteration: " << 0 << endl;
187 alpha[0] =
Blas::Ddot(ntot, &Kseqcopy[0][0], 1, &Kseqcopy[0][0], 1);
190 Vmath::Smul(ntot, 1.0/alpha[0], Kseq[0], 1, Kseq[0], 1);
194 for (i = 1; !converged && i <=
m_kdim; ++i)
204 alpha[i] =
Blas::Ddot(ntot, &Kseqcopy[i][0], 1, &Kseqcopy[i][0], 1);
209 Vmath::Smul(ntot, 1.0/alpha[i], Kseq[i], 1, Kseq[i], 1);
212 for (
int k = 0; k < i + 1; ++k)
226 EV_small(Tseq, Kseqcopy, ntot, alpha, i, zvec, wr, wi, resnorm);
229 converged =
EV_test(i, i, zvec, wr, wi, resnorm,
230 std::min(i,
m_nvec), evlout, resid0);
233 converged = max (converged, 0);
240 if (
m_comm->GetRank() == 0)
242 out <<
"Iteration: " << i <<
" (residual : " << resid0
254 for (
int j = 1; j <=
m_kdim; ++j)
256 alpha[j-1] = alpha[j];
277 for (
int k = 0; k <
m_kdim + 1; ++k)
291 EV_small(Tseq, Kseqcopy, ntot, alpha,
m_kdim, zvec, wr, wi, resnorm);
297 if (
m_comm->GetRank() == 0)
299 out <<
"Iteration: " << i <<
" (residual : "
300 << resid0 <<
")" <<endl;
311 for(j = 0; j <
m_equ[0]->GetNvariables(); ++j)
315 if (
m_comm->GetRank() == 0)
317 out <<
"L 2 error (variable " <<
m_equ[0]->GetVariable(j)
318 <<
") : " << vL2Error << endl;
319 out <<
"L inf error (variable " <<
m_equ[0]->GetVariable(j)
320 <<
") : " << vLinfError << endl;
330 "Need to implement Ritz re-evaluation of"
331 "eigenvalue. Only one half of complex "
332 "value will be correct");
339 if (
m_comm->GetRank() == 0)
355 m_equ[0]->TransCoeffToPhys();
357 m_equ[0]->SetTime(0.);
363 fields =
m_equ[0]->UpdateFields();
367 m_equ[1]->TransCoeffToPhys();
369 m_equ[1]->SetTime(0.);
393 int kdimp = kdim + 1;
401 for (
int i = 0; i < kdimp; ++i)
406 ASSERTL0(gsc != 0.0,
"Vectors are linearly independent.");
409 Vmath::Smul(ntot, 1.0/gsc, Kseq[i], 1, Kseq[i], 1);
412 Vmath::Smul(ntot, 1.0/gsc, Kseqcopy[i], 1, Kseqcopy[i], 1);
415 for (
int j = i + 1; j < kdimp; ++j)
417 gsc =
Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
419 Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
422 Vmath::Svtvp(ntot, -gsc, Kseqcopy[i], 1, Kseqcopy[j], 1, Kseqcopy[j], 1);
429 for (
int i = 0; i < kdim; ++i)
431 for (
int j = 0; j < kdim; ++j)
433 H[j*kdim+i] = alpha[j+1] * R[(j+1)*kdimp+i]
434 -
Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j*kdimp, 1);
435 H[j*kdim+i] /= R[j*kdimp+j];
439 H[(kdim-1)*kdim+kdim] = alpha[kdim]
440 * std::fabs(R[kdim*kdimp+kdim] / R[(kdim-1)*kdimp + kdim-1]);
442 Lapack::dgeev_(
'N',
'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1,
443 &zvec[0], kdim, &rwork[0], lwork, ier);
447 resnorm = H[(kdim-1)*kdim + kdim];
469 for (
int i = 0; i < kdim; ++i)
472 &zvec[0] + i*kdim, 1));
473 resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i*kdim]) / tmp;
476 resid[i-1] = resid[i] = hypot(resid[i-1], resid[i]);
480 EV_sort(zvec, wr, wi, resid, kdim);
482 while (nvec <= kdim && resid[nvec-1] <
m_evtol)
489 if (
m_comm->GetRank() == 0)
491 evlout <<
"-- Iteration = " << itrn <<
", H(k+1, k) = "
494 evlout.setf(ios::scientific, ios::floatfield);
497 evlout <<
" Magnitude Angle Growth "
498 <<
"Frequency Residual" << endl;
502 evlout <<
" Real Imaginary inverse real "
503 <<
"inverse imag Residual" << endl;
506 for (
int i = 0; i < kdim; i++)
508 WriteEvs(evlout,i,wr[i],wi[i],resid[i]);
512 resid0 = resid[nvec-1];
529 for (
int j = 1; j < dim; ++j)
536 while (i >= 0 && test[i] > te_tmp)
541 Vmath::Vcopy(dim, &evec[0] + i*dim, 1, &evec[0] + (i+1)*dim, 1);
547 Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i+1)*dim, 1);
571 ASSERTL0(
false,
"Convergence was not achieved within the "
572 "prescribed number of iterations.");
577 ASSERTL0(
false,
"Minimum residual reached.");
579 else if (icon == nvec)
582 EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
584 =
m_equ[0]->UpdateFields();
586 for (
int j = 0; j < icon; ++j)
588 std::string file =
m_session->GetSessionName() +
"_eig_"
589 + boost::lexical_cast<std::string>(j)
592 if (
m_comm->GetRank() == 0)
599 std::string fileunmask =
m_session->GetSessionName() +
"_eig_masked_"
600 + boost::lexical_cast<std::string>(j)
609 ASSERTL0(
false,
"Unrecognised value.");
627 boost::ignore_unused(wr);
632 for(
int i=0; i<nvec; ++i)
646 for (
int j = 0; j < nvec; ++j)
653 for (
int i = 0; i < kdim; ++i)
655 wgt = zvec[i + j*kdim];
656 Vmath::Svtvp(ntot, wgt, bvecs[i], 1, btmp[j], 1, btmp[j], 1);
659 Vmath::Svtvp(ntot, wgt, evecs[i], 1, etmp[j], 1, etmp[j], 1);
665 for (
int i = 0; i < nvec; ++i)
669 if (
m_session->GetComm()->GetRank() == 0)
671 cout <<
"eigenvalue " << i <<
": real mode" << endl;
673 norm =
Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
678 Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, bvecs[i], 1);
679 Vmath::Smul(ntot, 1.0/norm, etmp[i], 1, evecs[i], 1);
683 Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, evecs[i], 1);
688 if (
m_session->GetComm()->GetRank() == 0)
690 cout <<
"eigenvalues " << i <<
", " << i+1
691 <<
": complex modes" << endl;
693 norm =
Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
694 norm +=
Blas::Ddot(ntot, &btmp[i+1][0], 1, &btmp[i+1][0], 1);
700 Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, bvecs[i], 1);
701 Vmath::Smul(ntot, 1.0/norm, btmp[i+1], 1, bvecs[i+1], 1);
702 Vmath::Smul(ntot, 1.0/norm, etmp[i], 1, evecs[i], 1);
703 Vmath::Smul(ntot, 1.0/norm, etmp[i+1], 1, evecs[i+1], 1);
707 Vmath::Smul(ntot, 1.0/norm, btmp[i], 1, evecs[i], 1);
708 Vmath::Smul(ntot, 1.0/norm, btmp[i+1], 1, evecs[i+1], 1);
#define ASSERTL0(condition, msg)
#define WARNINGL0(condition, msg)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Base class for the development of solvers.
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)
Write coefficients to file.
SOLVER_UTILS_EXPORT const Array< OneD, const NekDouble > & GetMaskCoeff() const
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
virtual void v_InitObject(std::ostream &out=std::cout)
int m_nvec
Dimension of Krylov subspace.
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
int m_nits
Number of vectors to test.
Array< OneD, NekDouble > m_imag_evl
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
NekDouble m_evtol
Maxmum number of iterations.
int m_nfields
interval to dump information if required.
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
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)
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
LibUtilities::CommSharedPtr m_comm
Communication object.
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
int m_nequ
number of equations
virtual void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation.
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_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)
int EV_test(const int itrn, const int kdim, Array< OneD, NekDouble > &zvec, Array< OneD, NekDouble > &wr, Array< OneD, NekDouble > &wi, const NekDouble resnorm, int nvec, std::ofstream &evlout, NekDouble &resid0)
Tests for convergence of eigenvalues of H.
void EV_small(Array< OneD, Array< OneD, NekDouble > > &Kseq, Array< OneD, Array< OneD, NekDouble > > &Kseqcopy, 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 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_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)
virtual ~DriverModifiedArnoldi()
Destructor.
virtual void v_Execute(std::ostream &out=std::cout)
Virtual function for solve implementation.
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DriverFactory & GetDriverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The above copyright notice and this permission notice shall be included.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
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
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > sqrt(scalarT< T > in)