Nektar++
DriverArnoldi.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DriverArnoldi.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Base Driver class for the stability solver
32//
33///////////////////////////////////////////////////////////////////////////////
34#include <iomanip>
35
37
38using namespace std;
39
40namespace Nektar::SolverUtils
41{
42
43/**
44 * Constructor
45 */
49 : Driver(pSession, pGraph)
50{
51 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 1);
52}
53
54/**
55 * Arnoldi driver initialisation
56 */
58{
60 m_session->MatchSolverInfo("SolverType", "VelocityCorrectionScheme",
62
64 {
65 m_period = m_session->GetParameter("TimeStep") *
66 m_session->GetParameter("NumSteps");
67 m_nfields = m_equ[0]->UpdateFields().size() - 1;
68 }
69 else
70 {
71 m_period = 1.0;
72 m_nfields = m_equ[0]->UpdateFields().size();
73 }
74
75 if (m_session->DefinesSolverInfo("ModeType") &&
76 (boost::iequals(m_session->GetSolverInfo("ModeType"), "SingleMode") ||
77 boost::iequals(m_session->GetSolverInfo("ModeType"), "HalfMode")))
78 {
79 for (int i = 0; i < m_nfields; ++i)
80 {
81 m_equ[0]->UpdateFields()[i]->SetWaveSpace(true);
82 }
83 }
84 m_negatedOp = m_equ[0]->NegatedOp();
85
86 m_session->LoadParameter("kdim", m_kdim, 16);
87 m_session->LoadParameter("nvec", m_nvec, 2);
88 m_session->LoadParameter("nits", m_nits, 500);
89 m_session->LoadParameter("evtol", m_evtol, 1e-06);
90
91 ASSERTL0(m_kdim >= m_nvec, "nvec cannot be larger than kdim.");
92
93 m_session->LoadParameter("realShift", m_realShift, 0.0);
94 m_equ[0]->SetLambda(m_realShift);
95
96 m_session->LoadParameter("imagShift", m_imagShift, 0.0);
97
98 // The imaginary shift is applied at system assembly
99 // Only if using HOMOGENEOUS expansion and ModeType set to SingleMode
100 if (m_imagShift != 0.0)
101 {
102 if (!m_session->DefinesSolverInfo("HOMOGENEOUS") &&
103 !m_session->DefinesSolverInfo("ModeType"))
104 {
106 "Imaginary shift only supported with HOMOGENEOUS "
107 "expansion and ModeType set to SingleMode");
108 }
109 else if (!boost::iequals(m_session->GetSolverInfo("ModeType"),
110 "SingleMode"))
111 {
113 "Imaginary shift only supported with HOMOGENEOUS "
114 "expansion and ModeType set to SingleMode");
115 }
116 }
117 MaskInit();
118}
119
120/**
121 *
122 */
123void DriverArnoldi::v_Execute([[maybe_unused]] std::ostream &out)
124{
125 ASSERTL0(false, "Specific version of Arnoldi driver not implemented");
126}
127
128/**
129 *
130 */
131void DriverArnoldi::ArnoldiSummary(std::ostream &out)
132{
133 if (m_comm->GetRank() == 0)
134 {
135 if (m_session->DefinesSolverInfo("ModeType") &&
136 boost::iequals(m_session->GetSolverInfo("ModeType"), "SingleMode"))
137 {
138 out << "\tSingle Fourier mode : true " << endl;
139 ASSERTL0(m_session->DefinesSolverInfo("Homogeneous"),
140 "Expected a homogeneous expansion to be defined "
141 "with single mode");
142 }
143 else
144 {
145 out << "\tSingle Fourier mode : false " << endl;
146 }
147 if (m_session->DefinesSolverInfo("BetaZero"))
148 {
149 out << "\tBeta set to Zero : true (overrides LHom)" << endl;
150 }
151 else
152 {
153 out << "\tBeta set to Zero : false " << endl;
154 }
155
157 {
158 out << "\tEvolution operator : "
159 << m_session->GetSolverInfo("EvolutionOperator") << endl;
160 }
161 else
162 {
163 out << "\tShift (Real,Imag) : " << m_realShift << ","
164 << m_imagShift << endl;
165 }
166 out << "\tKrylov-space dimension : " << m_kdim << endl;
167 out << "\tNumber of vectors : " << m_nvec << endl;
168 out << "\tMax iterations : " << m_nits << endl;
169 out << "\tEigenvalue tolerance : " << m_evtol << endl;
170 out << "======================================================" << endl;
171 }
172}
173
174/**
175 * Copy Arnoldi array to field variables which depend from
176 * either the m_fields or m_forces
177 */
179{
180
182 m_equ[0]->UpdateFields();
183 int nq = fields[0]->GetNcoeffs();
184
185 for (int k = 0; k < m_nfields; ++k)
186 {
187 Vmath::Vcopy(nq, &array[k * nq], 1, &fields[k]->UpdateCoeffs()[0], 1);
188 fields[k]->SetPhysState(false);
189 }
190}
191
192/**
193 * Copy field variables which depend from either the m_fields
194 * or m_forces array the Arnoldi array
195 */
197{
198
200
202 {
203 // This matters for the Adaptive SFD method because
204 // m_equ[1] is the nonlinear problem with non
205 // homogeneous BCs.
206 fields = m_equ[0]->UpdateFields();
207 }
208 else
209 {
210 fields = m_equ[m_nequ - 1]->UpdateFields();
211 }
212
213 for (int k = 0; k < m_nfields; ++k)
214 {
215 int nq = fields[0]->GetNcoeffs();
216 Vmath::Vcopy(nq, &fields[k]->GetCoeffs()[0], 1, &array[k * nq], 1);
217 fields[k]->SetPhysState(false);
218 }
219}
220
221/**
222 * Initialisation for the transient growth
223 */
225{
227
229 "Transient Growth non available for Coupled Solver");
230
231 fields = m_equ[0]->UpdateFields();
232 int nq = fields[0]->GetNcoeffs();
233
234 for (int k = 0; k < m_nfields; ++k)
235 {
236 Vmath::Vcopy(nq, &fields[k]->GetCoeffs()[0], 1,
237 &m_equ[1]->UpdateFields()[k]->UpdateCoeffs()[0], 1);
238 }
239}
240
241/**
242 *
243 */
244void DriverArnoldi::WriteFld(std::string file,
245 std::vector<Array<OneD, NekDouble>> coeffs)
246{
247
248 std::vector<std::string> variables(m_nfields);
249
250 ASSERTL1(coeffs.size() >= m_nfields, "coeffs is not of the correct length");
251 for (int i = 0; i < m_nfields; ++i)
252 {
253 variables[i] = m_equ[0]->GetVariable(i);
254 }
255
256 m_equ[0]->WriteFld(file, m_equ[0]->UpdateFields()[0], coeffs, variables);
257}
258
259/**
260 *
261 */
263{
264
265 std::vector<std::string> variables(m_nfields);
266 std::vector<Array<OneD, NekDouble>> fieldcoeffs(m_nfields);
267
268 int ncoeffs = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
269 ASSERTL1(coeffs.size() >= ncoeffs * m_nfields,
270 "coeffs is not of sufficient size");
271
272 for (int i = 0; i < m_nfields; ++i)
273 {
274 variables[i] = m_equ[0]->GetVariable(i);
275 fieldcoeffs[i] = coeffs + i * ncoeffs;
276 }
277
278 m_equ[0]->WriteFld(file, m_equ[0]->UpdateFields()[0], fieldcoeffs,
279 variables);
280}
281
282/**
283 *
284 */
285void DriverArnoldi::WriteEvs(ostream &evlout, const int i,
286 const NekDouble re_ev, const NekDouble im_ev,
287 NekDouble resid, bool DumpInverse)
288{
290 {
291 NekDouble abs_ev = hypot(re_ev, im_ev);
292 NekDouble ang_ev = atan2(im_ev, re_ev);
293
294 evlout << "EV: " << setw(2) << i << setw(12) << abs_ev << setw(12)
295 << ang_ev << setw(12) << log(abs_ev) / m_period << setw(12)
296 << ang_ev / m_period;
297
299 {
300 evlout << setw(12) << resid;
301 }
302 evlout << endl;
303 }
304 else
305 {
306 NekDouble invmag = 1.0 / (re_ev * re_ev + im_ev * im_ev);
308 if (m_negatedOp)
309 {
310 sign = -1.0;
311 }
312 else
313 {
314 sign = 1.0;
315 }
316
317 evlout << "EV: " << setw(2) << i << setw(14) << sign * re_ev << setw(14)
318 << sign * im_ev;
319
320 if (DumpInverse)
321 {
322 evlout << setw(14) << sign * re_ev * invmag + m_realShift
323 << setw(14) << sign * im_ev * invmag + m_imagShift;
324 }
325
327 {
328 evlout << setw(12) << resid;
329 }
330 evlout << endl;
331 }
332}
333
334/**
335 *
336 */
338 std::vector<std::vector<LibUtilities::EquationSharedPtr>> &selectedDomains,
339 std::set<int> &unselectedVariables)
340{
341 selectedDomains.clear();
342 string domain("SelectEVCalcDomain0");
343 string condition("C0");
344 for (size_t i = 0; i < 10; ++i)
345 {
346 domain[domain.size() - 1] = '0' + i;
347 if (!m_session->DefinesFunction(domain))
348 {
349 break;
350 }
351 for (size_t j = 0; j < 10; ++j)
352 {
353 condition[condition.size() - 1] = '0' + j;
354 if (!m_session->DefinesFunction(domain, condition))
355 {
356 break;
357 }
358 if (j == 0)
359 {
360 selectedDomains.push_back(
361 std::vector<LibUtilities::EquationSharedPtr>());
362 }
363 selectedDomains[selectedDomains.size() - 1].push_back(
364 m_session->GetFunction(domain, condition));
365 }
366 }
367 unselectedVariables.clear();
368 string funName("SelectEVCalcVariables");
369 std::vector<std::string> variables = m_session->GetVariables();
370 for (size_t v = 0; v < m_nfields; ++v)
371 {
372 if (!m_session->DefinesFunction(funName, variables[v]))
373 {
374 unselectedVariables.insert(v);
375 }
376 }
377 if (unselectedVariables.size() == m_nfields)
378 {
379 unselectedVariables.clear();
380 }
381}
382
383/**
384 *
385 */
387{
388 std::vector<std::vector<LibUtilities::EquationSharedPtr>> selectedDomains;
389 std::set<int> unselectedVariables;
390 GetMaskInfo(selectedDomains, unselectedVariables);
391 if (selectedDomains.size() == 0 && unselectedVariables.size() == 0)
392 {
393 m_useMask = false;
394 return;
395 }
396 m_useMask = true;
397 MultiRegions::ExpListSharedPtr field = m_equ[0]->UpdateFields()[0];
398 int ncoef = field->GetNcoeffs();
399 int nphys = field->GetNpoints();
402 for (size_t i = 0; i < field->GetExpSize(); ++i)
403 {
404 LocalRegions::ExpansionSharedPtr exp = field->GetExp(i);
405 SpatialDomains::GeometrySharedPtr geom = exp->GetGeom();
406 int nv = geom->GetNumVerts();
407 NekDouble gc[3] = {0., 0., 0.};
408 NekDouble gct[3] = {0., 0., 0.};
409 for (size_t j = 0; j < nv; ++j)
410 {
411 SpatialDomains::PointGeomSharedPtr vertex = geom->GetVertex(j);
412 vertex->GetCoords(gct[0], gct[1], gct[2]);
413 gc[0] += gct[0] / NekDouble(nv);
414 gc[1] += gct[1] / NekDouble(nv);
415 gc[2] += gct[2] / NekDouble(nv);
416 }
417 int unmask = 1;
418 for (size_t m = 0; m < selectedDomains.size(); ++m)
419 {
420 unmask = 0;
421 for (size_t n = 0; n < selectedDomains[m].size(); ++n)
422 {
423 if (selectedDomains[m][n]->Evaluate(gc[0], gc[1], gc[2]) <= 0.)
424 {
425 unmask = 0;
426 break;
427 }
428 else
429 {
430 unmask = 1;
431 }
432 }
433 if (unmask == 1)
434 {
435 break;
436 }
437 }
438 for (int j = 0; j < m_nfields; ++j)
439 {
440 if (unmask == 0 || unselectedVariables.count(j))
441 {
443 exp->GetNcoeffs(), 0.,
444 &m_maskCoeffs[field->GetCoeff_Offset(i) + j * ncoef], 1);
445 Vmath::Fill(exp->GetTotPoints(), 0.,
446 &m_maskPhys[field->GetPhys_Offset(i) + j * nphys],
447 1);
448 }
449 }
450 }
451}
452
453} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
void v_InitObject(std::ostream &out=std::cout) override
Virtual function for initialisation implementation.
void WriteFld(std::string file, std::vector< Array< OneD, NekDouble > > coeffs)
Write coefficients to file.
NekDouble m_period
Tolerance of iterations.
Definition: DriverArnoldi.h:64
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:67
void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
Array< OneD, NekDouble > m_maskCoeffs
Definition: DriverArnoldi.h:78
DriverArnoldi(const LibUtilities::SessionReaderSharedPtr pSession, const SpatialDomains::MeshGraphSharedPtr pGraph)
Constructor.
void CopyFieldToArnoldiArray(Array< OneD, NekDouble > &array)
Copy fields to Arnoldi storage.
int m_nvec
Dimension of Krylov subspace.
Definition: DriverArnoldi.h:61
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:65
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:62
void CopyArnoldiArrayToField(Array< OneD, NekDouble > &array)
Copy Arnoldi storage to fields.
void GetMaskInfo(std::vector< std::vector< LibUtilities::EquationSharedPtr > > &selectedDomains, std::set< int > &unselectedVariables)
Array< OneD, NekDouble > m_maskPhys
Definition: DriverArnoldi.h:79
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:63
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:69
SOLVER_UTILS_EXPORT void ArnoldiSummary(std::ostream &out)
void WriteEvs(std::ostream &evlout, const int k, const NekDouble real, const NekDouble imag, NekDouble resid=NekConstants::kNekUnsetDouble, bool DumpInverse=true)
Base class for the development of solvers.
Definition: Driver.h:65
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:83
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation.
Definition: Driver.cpp:91
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:80
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:98
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:92
int m_nequ
number of equations
Definition: Driver.h:95
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kNekUnsetDouble
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:51
double NekDouble
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303