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...
 
virtual ~DriverModifiedArnoldi ()=default
 Destructor. More...
 
virtual void v_InitObject (std::ostream &out=std::cout) override
 Virtual function for initialisation implementation. More...
 
virtual 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...
 
virtual ~DriverArnoldi ()=default
 Destructor. More...
 
virtual void v_InitObject (std::ostream &out=std::cout) override
 Virtual function for initialisation implementation. More...
 
virtual 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 46 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 59 of file DriverModifiedArnoldi.cpp.

62 : DriverArnoldi(pSession, pGraph)
63{
64}
DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.

◆ ~DriverModifiedArnoldi()

virtual Nektar::SolverUtils::DriverModifiedArnoldi::~DriverModifiedArnoldi ( )
protectedvirtualdefault

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 52 of file DriverModifiedArnoldi.h.

55 {
58 pGraph);
59 p->InitObject();
60 return p;
61 }
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:54

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

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

544{
545 if (icon == 0)
546 {
547 // Not converged, write final Krylov vector
548 ASSERTL0(false, "Convergence was not achieved within the "
549 "prescribed number of iterations.");
550 }
551 else if (icon < 0)
552 {
553 // Minimum residual reached
554 ASSERTL0(false, "Minimum residual reached.");
555 }
556 else if (icon == nvec)
557 {
558 // Converged, write out eigenvectors
559 EV_big(Tseq, Kseq, ntot, kdim, icon, zvec, wr, wi);
560 Array<OneD, MultiRegions::ExpListSharedPtr> fields =
561 m_equ[0]->UpdateFields();
562
563 for (int j = 0; j < icon; ++j)
564 {
565 std::string file = m_session->GetSessionName() + "_eig_" +
566 boost::lexical_cast<std::string>(j) + ".fld";
567
568 if (m_comm->GetRank() == 0)
569 {
570 WriteEvs(cout, j, wr[j], wi[j]);
571 }
572 WriteFld(file, Kseq[j]);
573 if (m_useMask)
574 {
575 std::string fileunmask =
576 m_session->GetSessionName() + "_eig_masked_" +
577 boost::lexical_cast<std::string>(j) + ".fld";
578 WriteFld(fileunmask, Tseq[j]);
579 }
580 }
581 }
582 else
583 {
584 // Not recognised
585 ASSERTL0(false, "Unrecognised value.");
586 }
587}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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:94
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 371 of file DriverModifiedArnoldi.cpp.

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

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

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

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

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

342{
343 // Copy starting vector into first sequence element.
345 m_equ[0]->TransCoeffToPhys();
346
347 m_equ[0]->SetTime(0.);
348 m_equ[0]->DoSolve();
349
351 {
352 Array<OneD, MultiRegions::ExpListSharedPtr> fields;
353 fields = m_equ[0]->UpdateFields();
354
355 // start Adjoint with latest fields of direct
356 CopyFwdToAdj();
357 m_equ[1]->TransCoeffToPhys();
358
359 m_equ[1]->SetTime(0.);
360 m_equ[1]->DoSolve();
361 }
362
363 // Copy starting vector into first sequence element.
365}
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:100

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

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

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

70{
72
73 m_equ[0]->PrintSummary(out);
74
75 // Print session parameters
76 if (m_comm->GetRank() == 0)
77 {
78 out << "\tArnoldi solver type : Modified Arnoldi" << endl;
79 }
80
82
83 for (int i = 0; i < m_nequ; ++i)
84 {
85 m_equ[i]->DoInitialise();
86 }
87
88 // FwdTrans Initial conditions to be in Coefficient Space
89 m_equ[m_nequ - 1]->TransPhysToCoeff();
90}
virtual 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:97

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:198
static DriverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:67

Name of the class.

Definition at line 64 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 80 of file DriverModifiedArnoldi.h.