Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Static Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::SolverUtils::DriverModifiedArnoldi Class Reference

#include <DriverModifiedArnoldi.h>

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

Static Public Member Functions

static DriverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

 DriverModifiedArnoldi (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Constructor. More...
 
 ~DriverModifiedArnoldi () override=default
 Destructor. More...
 
void v_InitObject (std::ostream &out=std::cout) override
 Virtual function for initialisation implementation. More...
 
void v_Execute (std::ostream &out=std::cout) override
 Virtual function for solve implementation. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
 DriverArnoldi (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Constructor. More...
 
 ~DriverArnoldi () override=default
 Destructor. More...
 
void v_InitObject (std::ostream &out=std::cout) override
 Virtual function for initialisation implementation. More...
 
void v_Execute (std::ostream &out=std::cout) override
 Virtual function for solve implementation. More...
 
void CopyArnoldiArrayToField (Array< OneD, NekDouble > &array)
 Copy Arnoldi storage to fields. More...
 
void CopyFieldToArnoldiArray (Array< OneD, NekDouble > &array)
 Copy fields to Arnoldi storage. More...
 
void CopyFwdToAdj ()
 Copy the forward field to the adjoint system in transient growth calculations. More...
 
void WriteFld (std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
 Write coefficients to file. More...
 
void WriteFld (std::string file, Array< OneD, NekDouble > coeffs)
 
void WriteEvs (std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
 
void MaskInit ()
 Init mask. More...
 
void GetMaskInfo (std::vector< std::vector< LibUtilities::EquationSharedPtr > > &selectedDomains, std::set< int > &unselectedVariables)
 
SOLVER_UTILS_EXPORT void ArnoldiSummary (std::ostream &out)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Driver
 Driver (const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (std::ostream &out=std::cout)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_Execute (std::ostream &out=std::cout)=0
 Virtual function for solve implementation. More...
 

Static Protected Attributes

static std::string driverLookupId
 
- Static Protected Attributes inherited from Nektar::SolverUtils::Driver
static std::string evolutionOperatorLookupIds []
 
static std::string evolutionOperatorDef
 
static std::string driverDefault
 

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. More...
 
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. More...
 
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. More...
 
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. More...
 
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)
 

Friends

class MemoryManager< DriverModifiedArnoldi >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::DriverArnoldi
Array< OneD, NekDoubleGetRealEvl (void)
 
Array< OneD, NekDoubleGetImagEvl (void)
 
- Public Member Functions inherited from Nektar::SolverUtils::Driver
virtual ~Driver ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (std::ostream &out=std::cout)
 Initialise Object. More...
 
SOLVER_UTILS_EXPORT void Execute (std::ostream &out=std::cout)
 Execute driver. More...
 
SOLVER_UTILS_EXPORT Array< OneD, EquationSystemSharedPtrGetEqu ()
 
- Protected Attributes inherited from Nektar::SolverUtils::DriverArnoldi
int m_kdim
 
int m_nvec
 Dimension of Krylov subspace. More...
 
int m_nits
 Number of vectors to test. More...
 
NekDouble m_evtol
 Maxmum number of iterations. More...
 
NekDouble m_period
 Tolerance of iterations. More...
 
bool m_timeSteppingAlgorithm
 Period of time stepping algorithm. More...
 
int m_infosteps
 underlying operator is time stepping More...
 
int m_nfields
 interval to dump information if required. More...
 
NekDouble m_realShift
 
NekDouble m_imagShift
 
int m_negatedOp
 
Array< OneD, NekDoublem_real_evl
 Operator in solve call is negated. More...
 
Array< OneD, NekDoublem_imag_evl
 
bool m_useMask
 
Array< OneD, NekDoublem_maskCoeffs
 
Array< OneD, NekDoublem_maskPhys
 
- Protected Attributes inherited from Nektar::SolverUtils::Driver
LibUtilities::CommSharedPtr m_comm
 Communication object. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 Session reader object. More...
 
LibUtilities::SessionReaderSharedPtr session_LinNS
 Coupling between SFD and arnoldi. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 MeshGraph object. More...
 
Array< OneD, EquationSystemSharedPtrm_equ
 Equation system to solve. More...
 
int m_nequ
 number of equations More...
 
enum EvolutionOperatorType m_EvolutionOperator
 Evolution Operator. More...
 

Detailed Description

Definition at line 44 of file DriverModifiedArnoldi.h.

Constructor & Destructor Documentation

◆ DriverModifiedArnoldi()

Nektar::SolverUtils::DriverModifiedArnoldi::DriverModifiedArnoldi ( const LibUtilities::SessionReaderSharedPtr  pSession,
const SpatialDomains::MeshGraphSharedPtr  pGraph 
)
protected

Constructor.

Definition at line 55 of file DriverModifiedArnoldi.cpp.

58 : DriverArnoldi(pSession, pGraph)
59{
60}
DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.

◆ ~DriverModifiedArnoldi()

Nektar::SolverUtils::DriverModifiedArnoldi::~DriverModifiedArnoldi ( )
overrideprotecteddefault

Destructor.

Member Function Documentation

◆ create()

static DriverSharedPtr Nektar::SolverUtils::DriverModifiedArnoldi::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file DriverModifiedArnoldi.h.

53 {
56 pGraph);
57 p->InitObject();
58 return p;
59 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:52

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ EV_big()

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

Compute estimates of the eigenvalues/eigenvectors of the original system.

Definition at line 588 of file DriverModifiedArnoldi.cpp.

594{
595 NekDouble wgt, norm;
596 Array<OneD, Array<OneD, NekDouble>> btmp(nvec);
597 Array<OneD, Array<OneD, NekDouble>> etmp(nvec);
598 for (int i = 0; i < nvec; ++i)
599 {
600 if (m_useMask)
601 {
602 btmp[i] = Array<OneD, NekDouble>(ntot, 0.);
603 etmp[i] = Array<OneD, NekDouble>(ntot, 0.);
604 }
605 else
606 {
607 btmp[i] = evecs[i];
608 }
609 }
610
611 // Generate big eigenvectors
612 for (int j = 0; j < nvec; ++j)
613 {
614 Vmath::Zero(ntot, btmp[j], 1);
615 if (m_useMask)
616 {
617 Vmath::Zero(ntot, etmp[j], 1);
618 }
619 for (int i = 0; i < kdim; ++i)
620 {
621 wgt = zvec[i + j * kdim];
622 Vmath::Svtvp(ntot, wgt, bvecs[i], 1, btmp[j], 1, btmp[j], 1);
623 if (m_useMask)
624 {
625 Vmath::Svtvp(ntot, wgt, evecs[i], 1, etmp[j], 1, etmp[j], 1);
626 }
627 }
628 }
629
630 // Normalise the big eigenvectors
631 for (int i = 0; i < nvec; ++i)
632 {
633 if (wi[i] == 0.0) // Real mode
634 {
635 if (m_session->GetComm()->GetRank() == 0)
636 {
637 cout << "eigenvalue " << i << ": real mode" << endl;
638 }
639 norm = Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
640 m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
641 norm = std::sqrt(norm);
642 if (m_useMask)
643 {
644 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, bvecs[i], 1);
645 Vmath::Smul(ntot, 1.0 / norm, etmp[i], 1, evecs[i], 1);
646 }
647 else
648 {
649 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, evecs[i], 1);
650 }
651 }
652 else
653 {
654 if (m_session->GetComm()->GetRank() == 0)
655 {
656 cout << "eigenvalues " << i << ", " << i + 1
657 << ": complex modes" << endl;
658 }
659 norm = Blas::Ddot(ntot, &btmp[i][0], 1, &btmp[i][0], 1);
660 norm += Blas::Ddot(ntot, &btmp[i + 1][0], 1, &btmp[i + 1][0], 1);
661 m_comm->AllReduce(norm, Nektar::LibUtilities::ReduceSum);
662 norm = std::sqrt(norm);
663
664 if (m_useMask)
665 {
666 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, bvecs[i], 1);
667 Vmath::Smul(ntot, 1.0 / norm, btmp[i + 1], 1, bvecs[i + 1], 1);
668 Vmath::Smul(ntot, 1.0 / norm, etmp[i], 1, evecs[i], 1);
669 Vmath::Smul(ntot, 1.0 / norm, etmp[i + 1], 1, evecs[i + 1], 1);
670 }
671 else
672 {
673 Vmath::Smul(ntot, 1.0 / norm, btmp[i], 1, evecs[i], 1);
674 Vmath::Smul(ntot, 1.0 / norm, btmp[i + 1], 1, evecs[i + 1], 1);
675 }
676
677 i++;
678 }
679 }
680}
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:83
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:80
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:163
double NekDouble
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.
Definition: Vmath.hpp:396
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Blas::Ddot(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::Driver::m_session, Nektar::SolverUtils::DriverArnoldi::m_useMask, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), tinysimd::sqrt(), Vmath::Svtvp(), and Vmath::Zero().

Referenced by EV_post().

◆ EV_post()

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

Post-process the Ritz eigenvalues/eigenvectors of the matrix H, to compute estimations of the leading eigenvalues and eigenvectors of the original matrix.

Definition at line 533 of file DriverModifiedArnoldi.cpp.

540{
541 if (icon == 0)
542 {
543 // Not converged, write final Krylov vector
544 ASSERTL0(false, "Convergence was not achieved within the "
545 "prescribed number of iterations.");
546 }
547 else if (icon < 0)
548 {
549 // Minimum residual reached
550 ASSERTL0(false, "Minimum residual reached.");
551 }
552 else if (icon == nvec)
553 {
554 // Converged, write out eigenvectors
555 EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
556 Array<OneD, MultiRegions::ExpListSharedPtr> fields =
557 m_equ[0]->UpdateFields();
558
559 for (int j = 0; j < icon; ++j)
560 {
561 std::string file = m_session->GetSessionName() + "_eig_" +
562 std::to_string(j) + ".fld";
563
564 if (m_comm->GetRank() == 0)
565 {
566 WriteEvs(cout, j, wr[j], wi[j]);
567 }
568 WriteFld(file, Kseq[j]);
569 if (m_useMask)
570 {
571 std::string fileunmask = m_session->GetSessionName() +
572 "_eig_masked_" + std::to_string(j) +
573 ".fld";
574 WriteFld(fileunmask, Tseq[j]);
575 }
576 }
577 }
578 else
579 {
580 // Not recognised
581 ASSERTL0(false, "Unrecognised value.");
582 }
583}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:92
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)

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

Referenced by v_Execute().

◆ EV_small()

void Nektar::SolverUtils::DriverModifiedArnoldi::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 
)
private

Generates the upper Hessenberg matrix H and computes its eigenvalues.

Computes the Ritz eigenvalues and eigenvectors of the Hessenberg matrix constructed from the Krylov sequence.

Definition at line 367 of file DriverModifiedArnoldi.cpp.

373{
374 int kdimp = kdim + 1;
375 int lwork = 10 * kdim;
376 int ier;
377 Array<OneD, NekDouble> R(kdimp * kdimp, 0.0);
378 Array<OneD, NekDouble> H(kdimp * kdim, 0.0);
379 Array<OneD, NekDouble> rwork(lwork, 0.0);
380
381 // Modified G-S orthonormalisation
382 for (int i = 0; i < kdimp; ++i)
383 {
384 NekDouble gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[i][0], 1);
386 gsc = std::sqrt(gsc);
387 ASSERTL0(gsc != 0.0, "Vectors are linearly independent.");
388
389 R[i * kdimp + i] = gsc;
390 Vmath::Smul(ntot, 1.0 / gsc, Kseq[i], 1, Kseq[i], 1);
391 if (m_useMask)
392 {
393 Vmath::Smul(ntot, 1.0 / gsc, Kseqcopy[i], 1, Kseqcopy[i], 1);
394 }
395
396 for (int j = i + 1; j < kdimp; ++j)
397 {
398 gsc = Blas::Ddot(ntot, &Kseq[i][0], 1, &Kseq[j][0], 1);
400 Vmath::Svtvp(ntot, -gsc, Kseq[i], 1, Kseq[j], 1, Kseq[j], 1);
401 if (m_useMask)
402 {
403 Vmath::Svtvp(ntot, -gsc, Kseqcopy[i], 1, Kseqcopy[j], 1,
404 Kseqcopy[j], 1);
405 }
406 R[j * kdimp + i] = gsc;
407 }
408 }
409
410 // Compute H matrix
411 for (int i = 0; i < kdim; ++i)
412 {
413 for (int j = 0; j < kdim; ++j)
414 {
415 H[j * kdim + i] =
416 alpha[j + 1] * R[(j + 1) * kdimp + i] -
417 Vmath::Dot(j, &H[0] + i, kdim, &R[0] + j * kdimp, 1);
418 H[j * kdim + i] /= R[j * kdimp + j];
419 }
420 }
421
422 H[(kdim - 1) * kdim + kdim] =
423 alpha[kdim] *
424 std::fabs(R[kdim * kdimp + kdim] / R[(kdim - 1) * kdimp + kdim - 1]);
425
426 Lapack::dgeev_('N', 'V', kdim, &H[0], kdim, &wr[0], &wi[0], nullptr, 1,
427 &zvec[0], kdim, &rwork[0], lwork, ier);
428
429 ASSERTL0(!ier, "Error with dgeev");
430
431 resnorm = H[(kdim - 1) * kdim + kdim];
432}
T Dot(int n, const T *w, const T *x)
dot product
Definition: Vmath.hpp:761

References ASSERTL0, Blas::Ddot(), Vmath::Dot(), Nektar::SolverUtils::Driver::m_comm, Nektar::SolverUtils::DriverArnoldi::m_useMask, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), tinysimd::sqrt(), and Vmath::Svtvp().

Referenced by v_Execute().

◆ EV_sort()

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.

Sorts the computed eigenvalues by smallest residual first.

Definition at line 498 of file DriverModifiedArnoldi.cpp.

502{
503 Array<OneD, NekDouble> z_tmp(dim, 0.0);
504 NekDouble wr_tmp, wi_tmp, te_tmp;
505 for (int j = 1; j < dim; ++j)
506 {
507 wr_tmp = wr[j];
508 wi_tmp = wi[j];
509 te_tmp = test[j];
510 Vmath::Vcopy(dim, &evec[0] + j * dim, 1, &z_tmp[0], 1);
511 int i = j - 1;
512 while (i >= 0 && test[i] > te_tmp)
513 {
514 wr[i + 1] = wr[i];
515 wi[i + 1] = wi[i];
516 test[i + 1] = test[i];
517 Vmath::Vcopy(dim, &evec[0] + i * dim, 1, &evec[0] + (i + 1) * dim,
518 1);
519 i--;
520 }
521 wr[i + 1] = wr_tmp;
522 wi[i + 1] = wi_tmp;
523 test[i + 1] = te_tmp;
524 Vmath::Vcopy(dim, &z_tmp[0], 1, &evec[0] + (i + 1) * dim, 1);
525 }
526}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::UnitTests::test(), and Vmath::Vcopy().

Referenced by EV_test().

◆ EV_test()

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,
int  nvec,
std::ofstream &  evlout,
NekDouble resid0 
)
private

Tests for convergence of eigenvalues of H.

Check for convergence of the residuals of the eigenvalues of H.

Definition at line 437 of file DriverModifiedArnoldi.cpp.

443{
444 int idone = 0;
445 // NekDouble period = 0.1;
446
447 Array<OneD, NekDouble> resid(kdim);
448 for (int i = 0; i < kdim; ++i)
449 {
450 NekDouble tmp = std::sqrt(
451 Vmath::Dot(kdim, &zvec[0] + i * kdim, 1, &zvec[0] + i * kdim, 1));
452 resid[i] = resnorm * std::fabs(zvec[kdim - 1 + i * kdim]) / tmp;
453 if (wi[i] < 0.0)
454 {
455 resid[i - 1] = resid[i] = hypot(resid[i - 1], resid[i]);
456 }
457 }
458
459 EV_sort(zvec, wr, wi, resid, kdim);
460
461 while (nvec <= kdim && resid[nvec - 1] < m_evtol)
462 {
463 idone = nvec;
464 ++nvec;
465 }
466 nvec -= (idone > 0);
467
468 if (m_comm->GetRank() == 0)
469 {
470 evlout << "-- Iteration = " << itrn << ", H(k+1, k) = " << resnorm
471 << endl;
472 evlout.precision(4);
473 evlout.setf(ios::scientific, ios::floatfield);
475 {
476 evlout << " Magnitude Angle Growth "
477 << "Frequency Residual" << endl;
478 }
479 else
480 {
481 evlout << " Real Imaginary inverse real "
482 << "inverse imag Residual" << endl;
483 }
484
485 for (int i = 0; i < kdim; i++)
486 {
487 WriteEvs(evlout, i, wr[i], wi[i], resid[i]);
488 }
489 }
490
491 resid0 = resid[nvec - 1];
492 return idone;
493}
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:65
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:63
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.

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

Referenced by v_Execute().

◆ EV_update()

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 336 of file DriverModifiedArnoldi.cpp.

338{
339 // Copy starting vector into first sequence element.
341 m_equ[0]->TransCoeffToPhys();
342
343 m_equ[0]->SetTime(0.);
344 m_equ[0]->DoSolve();
345
347 {
348 Array<OneD, MultiRegions::ExpListSharedPtr> fields;
349 fields = m_equ[0]->UpdateFields();
350
351 // start Adjoint with latest fields of direct
352 CopyFwdToAdj();
353 m_equ[1]->TransCoeffToPhys();
354
355 m_equ[1]->SetTime(0.);
356 m_equ[1]->DoSolve();
357 }
358
359 // Copy starting vector into first sequence element.
361}
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:98

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().

◆ v_Execute()

void Nektar::SolverUtils::DriverModifiedArnoldi::v_Execute ( std::ostream &  out = std::cout)
overrideprotectedvirtual

Virtual function for solve implementation.

Reimplemented from Nektar::SolverUtils::DriverArnoldi.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

Definition at line 91 of file DriverModifiedArnoldi.cpp.

92{
93 int i = 0;
94 int j = 0;
95 int nq = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
96 int ntot = m_nfields * nq;
97 int converged = 0;
98 NekDouble resnorm = 0.0;
99 ofstream evlout;
100 std::string evlFile = m_session->GetSessionName() + ".evl";
101
102 if (m_comm->GetRank() == 0)
103 {
104 evlout.open(evlFile.c_str());
105 }
106
107 // Allocate memory
108 Array<OneD, NekDouble> alpha = Array<OneD, NekDouble>(m_kdim + 1, 0.0);
109 Array<OneD, NekDouble> wr = Array<OneD, NekDouble>(m_kdim, 0.0);
110 Array<OneD, NekDouble> wi = Array<OneD, NekDouble>(m_kdim, 0.0);
111 Array<OneD, NekDouble> zvec = Array<OneD, NekDouble>(m_kdim * m_kdim, 0.0);
112
113 Array<OneD, Array<OneD, NekDouble>> Kseq =
114 Array<OneD, Array<OneD, NekDouble>>(m_kdim + 1);
115 Array<OneD, Array<OneD, NekDouble>> Kseqcopy =
116 Array<OneD, Array<OneD, NekDouble>>(m_kdim + 1);
117 Array<OneD, Array<OneD, NekDouble>> Tseq =
118 Array<OneD, Array<OneD, NekDouble>>(m_kdim + 1);
119 for (i = 0; i < m_kdim + 1; ++i)
120 {
121 Kseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
122 if (m_useMask)
123 {
124 Kseqcopy[i] = Array<OneD, NekDouble>(ntot, 0.0);
125 }
126 else
127 {
128 Kseqcopy[i] = Kseq[i];
129 }
130 Tseq[i] = Array<OneD, NekDouble>(ntot, 0.0);
131 }
132
133 // Copy starting vector into second sequence element (temporary).
134 if (m_session->DefinesFunction("InitialConditions"))
135 {
136 if (m_comm->GetRank() == 0)
137 {
138 out << "\tInital vector : specified in input file " << endl;
139 }
140 m_equ[0]->SetInitialConditions(0.0, false);
141
143 }
144 else
145 {
146 if (m_comm->GetRank() == 0)
147 {
148 out << "\tInital vector : random " << endl;
149 }
150
151 NekDouble eps = 0.0001;
152 Vmath::FillWhiteNoise(ntot, eps, &Kseq[1][0], 1);
153 if (m_useMask)
154 {
155 Vmath::Vmul(ntot, Kseq[1], 1, m_maskCoeffs, 1, Kseq[1], 1);
156 }
157 }
158
159 // Perform one iteration to enforce boundary conditions.
160 // Set this as the initial value in the sequence.
161 EV_update(Kseq[1], Kseq[0]);
162 if (m_comm->GetRank() == 0)
163 {
164 out << "Iteration: " << 0 << endl;
165 }
166
167 // Normalise first vector in sequence
168 if (m_useMask)
169 {
170 Vmath::Vmul(ntot, Kseq[0], 1, m_maskCoeffs, 1, Kseqcopy[0], 1);
171 }
172 alpha[0] = Blas::Ddot(ntot, &Kseqcopy[0][0], 1, &Kseqcopy[0][0], 1);
173 m_comm->AllReduce(alpha[0], Nektar::LibUtilities::ReduceSum);
174 alpha[0] = std::sqrt(alpha[0]);
175 Vmath::Smul(ntot, 1.0 / alpha[0], Kseq[0], 1, Kseq[0], 1);
176
177 // Fill initial krylov sequence
178 NekDouble resid0;
179 for (i = 1; !converged && i <= m_kdim; ++i)
180 {
181 // Compute next vector
182 EV_update(Kseq[i - 1], Kseq[i]);
183
184 // Normalise
185 if (m_useMask)
186 {
187 Vmath::Vmul(ntot, Kseq[i], 1, m_maskCoeffs, 1, Kseqcopy[i], 1);
188 }
189 alpha[i] = Blas::Ddot(ntot, &Kseqcopy[i][0], 1, &Kseqcopy[i][0], 1);
190 m_comm->AllReduce(alpha[i], Nektar::LibUtilities::ReduceSum);
191 alpha[i] = std::sqrt(alpha[i]);
192
193 // alpha[i] = std::sqrt(alpha[i]);
194 Vmath::Smul(ntot, 1.0 / alpha[i], Kseq[i], 1, Kseq[i], 1);
195
196 // Copy Krylov sequence into temporary storage
197 for (int k = 0; k < i + 1; ++k)
198 {
199 if (m_useMask)
200 {
201 Vmath::Vmul(ntot, Kseq[k], 1, m_maskCoeffs, 1, Tseq[k], 1);
202 Vmath::Vcopy(ntot, Kseq[k], 1, Kseqcopy[k], 1);
203 }
204 else
205 {
206 Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
207 ;
208 }
209 }
210
211 // Generate Hessenberg matrix and compute eigenvalues of it.
212 EV_small(Tseq, Kseqcopy, ntot, alpha, i, zvec, wr, wi, resnorm);
213
214 // Test for convergence.
215 converged = EV_test(i, i, zvec, wr, wi, resnorm, std::min(i, m_nvec),
216 evlout, resid0);
217 if (i >= m_nvec)
218 {
219 converged = max(converged, 0);
220 }
221 else
222 {
223 converged = 0;
224 }
225
226 if (m_comm->GetRank() == 0)
227 {
228 out << "Iteration: " << i << " (residual : " << resid0 << ")"
229 << endl;
230 }
231 }
232
233 // Continue with full sequence
234 if (!converged)
235 {
236 for (i = m_kdim + 1; !converged && i <= m_nits; ++i)
237 {
238 // Shift all the vectors in the sequence.
239 // First vector is removed.
240 for (int j = 1; j <= m_kdim; ++j)
241 {
242 alpha[j - 1] = alpha[j];
243 Vmath::Vcopy(ntot, Kseq[j], 1, Kseq[j - 1], 1);
244 }
245
246 // Compute next vector
247 EV_update(Kseq[m_kdim - 1], Kseq[m_kdim]);
248
249 // Compute new scale factor
250 if (m_useMask)
251 {
252 Vmath::Vmul(ntot, Kseq[m_kdim], 1, m_maskCoeffs, 1,
253 Kseqcopy[m_kdim], 1);
254 }
255 alpha[m_kdim] = Blas::Ddot(ntot, &Kseqcopy[m_kdim][0], 1,
256 &Kseqcopy[m_kdim][0], 1);
258 alpha[m_kdim] = std::sqrt(alpha[m_kdim]);
259 Vmath::Smul(ntot, 1.0 / alpha[m_kdim], Kseq[m_kdim], 1,
260 Kseq[m_kdim], 1);
261
262 // Copy Krylov sequence into temporary storage
263 for (int k = 0; k < m_kdim + 1; ++k)
264 {
265 if (m_useMask)
266 {
267 Vmath::Vmul(ntot, Kseq[k], 1, m_maskCoeffs, 1, Tseq[k], 1);
268 Vmath::Vcopy(ntot, Kseq[k], 1, Kseqcopy[k], 1);
269 }
270 else
271 {
272 Vmath::Vcopy(ntot, Kseq[k], 1, Tseq[k], 1);
273 ;
274 }
275 }
276
277 // Generate Hessenberg matrix and compute eigenvalues of it
278 EV_small(Tseq, Kseqcopy, ntot, alpha, m_kdim, zvec, wr, wi,
279 resnorm);
280
281 // Test for convergence.
282 converged = EV_test(i, m_kdim, zvec, wr, wi, resnorm, m_nvec,
283 evlout, resid0);
284
285 if (m_comm->GetRank() == 0)
286 {
287 out << "Iteration: " << i << " (residual : " << resid0 << ")"
288 << endl;
289 }
290 }
291 }
292
293 m_equ[0]->Output();
294
295 // Evaluate and output computation time and solution accuracy.
296 // The specific format of the error output is essential for the
297 // regression tests to work.
298 // Evaluate L2 Error
299 for (j = 0; j < m_equ[0]->GetNvariables(); ++j)
300 {
301 NekDouble vL2Error = m_equ[0]->L2Error(j, false);
302 NekDouble vLinfError = m_equ[0]->LinfError(j);
303 if (m_comm->GetRank() == 0)
304 {
305 out << "L 2 error (variable " << m_equ[0]->GetVariable(j)
306 << ") : " << vL2Error << endl;
307 out << "L inf error (variable " << m_equ[0]->GetVariable(j)
308 << ") : " << vLinfError << endl;
309 }
310 }
311
312 // Process eigenvectors and write out.
313 m_nvec = converged;
314 EV_post(Tseq, Kseqcopy, ntot, min(--i, m_kdim), m_nvec, zvec, wr, wi,
315 converged);
316
317 WARNINGL0(m_imagShift == 0, "Complex Shift applied. "
318 "Need to implement Ritz re-evaluation of"
319 "eigenvalue. Only one half of complex "
320 "value will be correct");
321
322 // store eigenvalues so they can be accessed from driver class
323 m_real_evl = wr;
324 m_imag_evl = wi;
325
326 // Close the runtime info file.
327 if (m_comm->GetRank() == 0)
328 {
329 evlout.close();
330 }
331}
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
Array< OneD, NekDouble > m_maskCoeffs
Definition: DriverArnoldi.h:78
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:61
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:62
Array< OneD, NekDouble > m_imag_evl
Definition: DriverArnoldi.h:75
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:69
Array< OneD, NekDouble > m_real_evl
Operator in solve call is negated.
Definition: DriverArnoldi.h:74
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 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.
Definition: Vmath.hpp:72
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
Definition: Vmath.cpp:154

References Nektar::SolverUtils::DriverArnoldi::CopyFieldToArnoldiArray(), Blas::Ddot(), 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_imagShift, Nektar::SolverUtils::DriverArnoldi::m_kdim, Nektar::SolverUtils::DriverArnoldi::m_maskCoeffs, 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::SolverUtils::DriverArnoldi::m_useMask, Nektar::LibUtilities::ReduceSum, Vmath::Smul(), tinysimd::sqrt(), Vmath::Vcopy(), Vmath::Vmul(), and WARNINGL0.

Referenced by Nektar::SolverUtils::DriverSteadyState::v_Execute().

◆ v_InitObject()

void Nektar::SolverUtils::DriverModifiedArnoldi::v_InitObject ( std::ostream &  out = std::cout)
overrideprotectedvirtual

Virtual function for initialisation implementation.

Reimplemented from Nektar::SolverUtils::DriverArnoldi.

Reimplemented in Nektar::SolverUtils::DriverSteadyState.

Definition at line 65 of file DriverModifiedArnoldi.cpp.

66{
68
69 m_equ[0]->PrintSummary(out);
70
71 // Print session parameters
72 if (m_comm->GetRank() == 0)
73 {
74 out << "\tArnoldi solver type : Modified Arnoldi" << endl;
75 }
76
78
79 for (int i = 0; i < m_nequ; ++i)
80 {
81 m_equ[i]->DoInitialise();
82 }
83
84 // FwdTrans Initial conditions to be in Coefficient Space
85 m_equ[m_nequ - 1]->TransPhysToCoeff();
86}
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
int m_nequ
number of equations
Definition: Driver.h:95

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

Referenced by Nektar::SolverUtils::DriverSteadyState::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< DriverModifiedArnoldi >

friend class MemoryManager< DriverModifiedArnoldi >
friend

Definition at line 1 of file DriverModifiedArnoldi.h.

Member Data Documentation

◆ className

string Nektar::SolverUtils::DriverModifiedArnoldi::className
static
Initial value:
=
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:65

Name of the class.

Definition at line 62 of file DriverModifiedArnoldi.h.

◆ driverLookupId

string Nektar::SolverUtils::DriverModifiedArnoldi::driverLookupId
staticprotected
Initial value:
=
0)
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 78 of file DriverModifiedArnoldi.h.