38 #include <boost/core/ignore_unused.hpp>
49 string DriverModifiedArnoldi::className =
51 DriverModifiedArnoldi::create);
52 string DriverModifiedArnoldi::driverLookupId =
53 LibUtilities::SessionReader::RegisterEnumValue(
"Driver",
"ModifiedArnoldi",
59 DriverModifiedArnoldi::DriverModifiedArnoldi(
80 m_equ[0]->PrintSummary(out);
83 if (
m_comm->GetRank() == 0)
85 out <<
"\tArnoldi solver type : Modified Arnoldi" << endl;
90 for (
int i = 0; i <
m_nequ; ++i)
92 m_equ[i]->DoInitialise();
106 int nq =
m_equ[0]->UpdateFields()[0]->GetNcoeffs();
111 std::string evlFile =
m_session->GetSessionName() +
".evl";
113 if (
m_comm->GetRank() == 0)
115 evlout.open(evlFile.c_str());
130 for (i = 0; i <
m_kdim + 1; ++i)
139 Kseqcopy[i] = Kseq[i];
145 if (
m_session->DefinesFunction(
"InitialConditions"))
147 if (
m_comm->GetRank() == 0)
149 out <<
"\tInital vector : specified in input file " << endl;
151 m_equ[0]->SetInitialConditions(0.0,
false);
157 if (
m_comm->GetRank() == 0)
159 out <<
"\tInital vector : random " << endl;
173 if (
m_comm->GetRank() == 0)
175 out <<
"Iteration: " << 0 << endl;
183 alpha[0] =
Blas::Ddot(ntot, &Kseqcopy[0][0], 1, &Kseqcopy[0][0], 1);
186 Vmath::Smul(ntot, 1.0 / alpha[0], Kseq[0], 1, Kseq[0], 1);
190 for (i = 1; !converged && i <=
m_kdim; ++i)
200 alpha[i] =
Blas::Ddot(ntot, &Kseqcopy[i][0], 1, &Kseqcopy[i][0], 1);
205 Vmath::Smul(ntot, 1.0 / alpha[i], Kseq[i], 1, Kseq[i], 1);
208 for (
int k = 0; k < i + 1; ++k)
223 EV_small(Tseq, Kseqcopy, ntot, alpha, i, zvec, wr, wi, resnorm);
226 converged =
EV_test(i, i, zvec, wr, wi, resnorm, std::min(i,
m_nvec),
230 converged = max(converged, 0);
237 if (
m_comm->GetRank() == 0)
239 out <<
"Iteration: " << i <<
" (residual : " << resid0 <<
")"
251 for (
int j = 1; j <=
m_kdim; ++j)
253 alpha[j - 1] = alpha[j];
274 for (
int k = 0; k <
m_kdim + 1; ++k)
297 if (
m_comm->GetRank() == 0)
299 out <<
"Iteration: " << i <<
" (residual : " << resid0 <<
")"
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)
353 m_equ[0]->TransCoeffToPhys();
355 m_equ[0]->SetTime(0.);
361 fields =
m_equ[0]->UpdateFields();
365 m_equ[1]->TransCoeffToPhys();
367 m_equ[1]->SetTime(0.);
386 int kdimp = kdim + 1;
387 int lwork = 10 * kdim;
394 for (
int i = 0; i < kdimp; ++i)
399 ASSERTL0(gsc != 0.0,
"Vectors are linearly independent.");
401 R[i * kdimp + i] = gsc;
402 Vmath::Smul(ntot, 1.0 / gsc, Kseq[i], 1, Kseq[i], 1);
405 Vmath::Smul(ntot, 1.0 / gsc, Kseqcopy[i], 1, Kseqcopy[i], 1);
408 for (
int j = i + 1; j < kdimp; ++j)
410 gsc =
Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
412 Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
415 Vmath::Svtvp(ntot, -gsc, Kseqcopy[i], 1, Kseqcopy[j], 1,
418 R[j * kdimp + i] = gsc;
423 for (
int i = 0; i < kdim; ++i)
425 for (
int j = 0; j < kdim; ++j)
428 alpha[j + 1] * R[(j + 1) * kdimp + i] -
429 Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j * kdimp, 1);
430 H[j * kdim + i] /= R[j * kdimp + j];
434 H[(kdim - 1) * kdim + kdim] =
436 std::fabs(R[kdim * kdimp + kdim] / R[(kdim - 1) * kdimp + kdim - 1]);
438 Lapack::dgeev_(
'N',
'V', kdim, &H[0], kdim, &wr[0], &wi[0], 0, 1, &zvec[0],
439 kdim, &rwork[0], lwork, ier);
443 resnorm = H[(kdim - 1) * kdim + kdim];
460 for (
int i = 0; i < kdim; ++i)
463 Vmath::Dot(kdim, &zvec[0] + i * kdim, 1, &zvec[0] + i * kdim, 1));
464 resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i * kdim]) / tmp;
467 resid[i - 1] = resid[i] = hypot(resid[i - 1], resid[i]);
471 EV_sort(zvec, wr, wi, resid, kdim);
473 while (nvec <= kdim && resid[nvec - 1] <
m_evtol)
480 if (
m_comm->GetRank() == 0)
482 evlout <<
"-- Iteration = " << itrn <<
", H(k+1, k) = " << resnorm
485 evlout.setf(ios::scientific, ios::floatfield);
488 evlout <<
" Magnitude Angle Growth "
489 <<
"Frequency Residual" << endl;
493 evlout <<
" Real Imaginary inverse real "
494 <<
"inverse imag Residual" << endl;
497 for (
int i = 0; i < kdim; i++)
499 WriteEvs(evlout, i, wr[i], wi[i], resid[i]);
503 resid0 = resid[nvec - 1];
517 for (
int j = 1; j < dim; ++j)
524 while (i >= 0 && test[i] > te_tmp)
528 test[i + 1] = test[i];
529 Vmath::Vcopy(dim, &evec[0] + i * dim, 1, &evec[0] + (i + 1) * dim,
535 test[i + 1] = te_tmp;
536 Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i + 1) * dim, 1);
547 const int ntot,
const int kdim,
556 ASSERTL0(
false,
"Convergence was not achieved within the "
557 "prescribed number of iterations.");
562 ASSERTL0(
false,
"Minimum residual reached.");
564 else if (icon == nvec)
567 EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
569 m_equ[0]->UpdateFields();
571 for (
int j = 0; j < icon; ++j)
573 std::string file =
m_session->GetSessionName() +
"_eig_" +
574 boost::lexical_cast<std::string>(j) +
".fld";
576 if (
m_comm->GetRank() == 0)
583 std::string fileunmask =
584 m_session->GetSessionName() +
"_eig_masked_" +
585 boost::lexical_cast<std::string>(j) +
".fld";
593 ASSERTL0(
false,
"Unrecognised value.");
602 const int ntot,
const int kdim,
607 boost::ignore_unused(wr);
612 for (
int i = 0; i < nvec; ++i)
626 for (
int j = 0; j < nvec; ++j)
633 for (
int i = 0; i < kdim; ++i)
635 wgt = zvec[i + j * kdim];
636 Vmath::Svtvp(ntot, wgt, bvecs[i], 1, btmp[j], 1, btmp[j], 1);
639 Vmath::Svtvp(ntot, wgt, evecs[i], 1, etmp[j], 1, etmp[j], 1);
645 for (
int i = 0; i < nvec; ++i)
649 if (
m_session->GetComm()->GetRank() == 0)
651 cout <<
"eigenvalue " << i <<
": real mode" << endl;
653 norm =
Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
658 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, bvecs[i], 1);
659 Vmath::Smul(ntot, 1.0 / norm, etmp[i], 1, evecs[i], 1);
663 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, evecs[i], 1);
668 if (
m_session->GetComm()->GetRank() == 0)
670 cout <<
"eigenvalues " << i <<
", " << i + 1
671 <<
": complex modes" << endl;
673 norm =
Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
674 norm +=
Blas::Ddot(ntot, &btmp[i + 1][0], 1, &btmp[i + 1][0], 1);
680 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, bvecs[i], 1);
681 Vmath::Smul(ntot, 1.0 / norm, btmp[i + 1], 1, bvecs[i + 1], 1);
682 Vmath::Smul(ntot, 1.0 / norm, etmp[i], 1, evecs[i], 1);
683 Vmath::Smul(ntot, 1.0 / norm, etmp[i + 1], 1, evecs[i + 1], 1);
687 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, evecs[i], 1);
688 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.
virtual void v_InitObject(std::ostream &out=std::cout) override
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.
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
void EV_update(Array< OneD, NekDouble > &src, Array< OneD, NekDouble > &tgt)
Generates a new vector in the sequence by applying the linear operator.
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_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_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_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)
virtual void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
virtual void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
virtual ~DriverModifiedArnoldi()
Destructor.
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)
dot (vector times vector): z = w*x
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)