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
41{
42namespace SolverUtils
43{
44
45/**
46 * Constructor
47 */
51 : Driver(pSession, pGraph)
52{
53 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 1);
54}
55
56/**
57 * Arnoldi driver initialisation
58 */
60{
62 m_session->MatchSolverInfo("SolverType", "VelocityCorrectionScheme",
64
66 {
67 m_period = m_session->GetParameter("TimeStep") *
68 m_session->GetParameter("NumSteps");
69 m_nfields = m_equ[0]->UpdateFields().size() - 1;
70 }
71 else
72 {
73 m_period = 1.0;
74 m_nfields = m_equ[0]->UpdateFields().size();
75 }
76
77 if (m_session->DefinesSolverInfo("ModeType") &&
78 (boost::iequals(m_session->GetSolverInfo("ModeType"), "SingleMode") ||
79 boost::iequals(m_session->GetSolverInfo("ModeType"), "HalfMode")))
80 {
81 for (int i = 0; i < m_nfields; ++i)
82 {
83 m_equ[0]->UpdateFields()[i]->SetWaveSpace(true);
84 }
85 }
86 m_negatedOp = m_equ[0]->NegatedOp();
87
88 m_session->LoadParameter("kdim", m_kdim, 16);
89 m_session->LoadParameter("nvec", m_nvec, 2);
90 m_session->LoadParameter("nits", m_nits, 500);
91 m_session->LoadParameter("evtol", m_evtol, 1e-06);
92
93 ASSERTL0(m_kdim >= m_nvec, "nvec cannot be larger than kdim.");
94
95 m_session->LoadParameter("realShift", m_realShift, 0.0);
96 m_equ[0]->SetLambda(m_realShift);
97
98 m_session->LoadParameter("imagShift", m_imagShift, 0.0);
99
100 // The imaginary shift is applied at system assembly
101 // Only if using HOMOGENEOUS expansion and ModeType set to SingleMode
102 if (m_imagShift != 0.0)
103 {
104 if (!m_session->DefinesSolverInfo("HOMOGENEOUS") &&
105 !m_session->DefinesSolverInfo("ModeType"))
106 {
108 "Imaginary shift only supported with HOMOGENEOUS "
109 "expansion and ModeType set to SingleMode");
110 }
111 else if (!boost::iequals(m_session->GetSolverInfo("ModeType"),
112 "SingleMode"))
113 {
115 "Imaginary shift only supported with HOMOGENEOUS "
116 "expansion and ModeType set to SingleMode");
117 }
118 }
119 MaskInit();
120}
121
122/**
123 *
124 */
125void DriverArnoldi::v_Execute(std::ostream &out)
126{
127 boost::ignore_unused(out);
128
129 ASSERTL0(false, "Specific version of Arnoldi driver not implemented");
130}
131
132/**
133 *
134 */
135void DriverArnoldi::ArnoldiSummary(std::ostream &out)
136{
137 if (m_comm->GetRank() == 0)
138 {
139 if (m_session->DefinesSolverInfo("ModeType") &&
140 boost::iequals(m_session->GetSolverInfo("ModeType"), "SingleMode"))
141 {
142 out << "\tSingle Fourier mode : true " << endl;
143 ASSERTL0(m_session->DefinesSolverInfo("Homogeneous"),
144 "Expected a homogeneous expansion to be defined "
145 "with single mode");
146 }
147 else
148 {
149 out << "\tSingle Fourier mode : false " << endl;
150 }
151 if (m_session->DefinesSolverInfo("BetaZero"))
152 {
153 out << "\tBeta set to Zero : true (overrides LHom)" << endl;
154 }
155 else
156 {
157 out << "\tBeta set to Zero : false " << endl;
158 }
159
161 {
162 out << "\tEvolution operator : "
163 << m_session->GetSolverInfo("EvolutionOperator") << endl;
164 }
165 else
166 {
167 out << "\tShift (Real,Imag) : " << m_realShift << ","
168 << m_imagShift << endl;
169 }
170 out << "\tKrylov-space dimension : " << m_kdim << endl;
171 out << "\tNumber of vectors : " << m_nvec << endl;
172 out << "\tMax iterations : " << m_nits << endl;
173 out << "\tEigenvalue tolerance : " << m_evtol << endl;
174 out << "======================================================" << endl;
175 }
176}
177
178/**
179 * Copy Arnoldi array to field variables which depend from
180 * either the m_fields or m_forces
181 */
183{
184
186 m_equ[0]->UpdateFields();
187 int nq = fields[0]->GetNcoeffs();
188
189 for (int k = 0; k < m_nfields; ++k)
190 {
191 Vmath::Vcopy(nq, &array[k * nq], 1, &fields[k]->UpdateCoeffs()[0], 1);
192 fields[k]->SetPhysState(false);
193 }
194}
195
196/**
197 * Copy field variables which depend from either the m_fields
198 * or m_forces array the Arnoldi array
199 */
201{
202
204
206 {
207 // This matters for the Adaptive SFD method because
208 // m_equ[1] is the nonlinear problem with non
209 // homogeneous BCs.
210 fields = m_equ[0]->UpdateFields();
211 }
212 else
213 {
214 fields = m_equ[m_nequ - 1]->UpdateFields();
215 }
216
217 for (int k = 0; k < m_nfields; ++k)
218 {
219 int nq = fields[0]->GetNcoeffs();
220 Vmath::Vcopy(nq, &fields[k]->GetCoeffs()[0], 1, &array[k * nq], 1);
221 fields[k]->SetPhysState(false);
222 }
223}
224
225/**
226 * Initialisation for the transient growth
227 */
229{
231
233 "Transient Growth non available for Coupled Solver");
234
235 fields = m_equ[0]->UpdateFields();
236 int nq = fields[0]->GetNcoeffs();
237
238 for (int k = 0; k < m_nfields; ++k)
239 {
240 Vmath::Vcopy(nq, &fields[k]->GetCoeffs()[0], 1,
241 &m_equ[1]->UpdateFields()[k]->UpdateCoeffs()[0], 1);
242 }
243}
244
245/**
246 *
247 */
248void DriverArnoldi::WriteFld(std::string file,
249 std::vector<Array<OneD, NekDouble>> coeffs)
250{
251
252 std::vector<std::string> variables(m_nfields);
253
254 ASSERTL1(coeffs.size() >= m_nfields, "coeffs is not of the correct length");
255 for (int i = 0; i < m_nfields; ++i)
256 {
257 variables[i] = m_equ[0]->GetVariable(i);
258 }
259
260 m_equ[0]->WriteFld(file, m_equ[0]->UpdateFields()[0], coeffs, variables);
261}
262
263/**
264 *
265 */
267{
268
269 std::vector<std::string> variables(m_nfields);
270 std::vector<Array<OneD, NekDouble>> fieldcoeffs(m_nfields);
271
272 int ncoeffs = m_equ[0]->UpdateFields()[0]->GetNcoeffs();
273 ASSERTL1(coeffs.size() >= ncoeffs * m_nfields,
274 "coeffs is not of sufficient size");
275
276 for (int i = 0; i < m_nfields; ++i)
277 {
278 variables[i] = m_equ[0]->GetVariable(i);
279 fieldcoeffs[i] = coeffs + i * ncoeffs;
280 }
281
282 m_equ[0]->WriteFld(file, m_equ[0]->UpdateFields()[0], fieldcoeffs,
283 variables);
284}
285
286/**
287 *
288 */
289void DriverArnoldi::WriteEvs(ostream &evlout, const int i,
290 const NekDouble re_ev, const NekDouble im_ev,
291 NekDouble resid, bool DumpInverse)
292{
294 {
295 NekDouble abs_ev = hypot(re_ev, im_ev);
296 NekDouble ang_ev = atan2(im_ev, re_ev);
297
298 evlout << "EV: " << setw(2) << i << setw(12) << abs_ev << setw(12)
299 << ang_ev << setw(12) << log(abs_ev) / m_period << setw(12)
300 << ang_ev / m_period;
301
303 {
304 evlout << setw(12) << resid;
305 }
306 evlout << endl;
307 }
308 else
309 {
310 NekDouble invmag = 1.0 / (re_ev * re_ev + im_ev * im_ev);
312 if (m_negatedOp)
313 {
314 sign = -1.0;
315 }
316 else
317 {
318 sign = 1.0;
319 }
320
321 evlout << "EV: " << setw(2) << i << setw(14) << sign * re_ev << setw(14)
322 << sign * im_ev;
323
324 if (DumpInverse)
325 {
326 evlout << setw(14) << sign * re_ev * invmag + m_realShift
327 << setw(14) << sign * im_ev * invmag + m_imagShift;
328 }
329
331 {
332 evlout << setw(12) << resid;
333 }
334 evlout << endl;
335 }
336}
337
338/**
339 *
340 */
342 std::vector<std::vector<LibUtilities::EquationSharedPtr>> &selectedDomains,
343 std::set<int> &unselectedVariables)
344{
345 selectedDomains.clear();
346 string domain("SelectEVCalcDomain0");
347 string condition("C0");
348 for (size_t i = 0; i < 10; ++i)
349 {
350 domain[domain.size() - 1] = '0' + i;
351 if (!m_session->DefinesFunction(domain))
352 {
353 break;
354 }
355 for (size_t j = 0; j < 10; ++j)
356 {
357 condition[condition.size() - 1] = '0' + j;
358 if (!m_session->DefinesFunction(domain, condition))
359 {
360 break;
361 }
362 if (j == 0)
363 {
364 selectedDomains.push_back(
365 std::vector<LibUtilities::EquationSharedPtr>());
366 }
367 selectedDomains[selectedDomains.size() - 1].push_back(
368 m_session->GetFunction(domain, condition));
369 }
370 }
371 unselectedVariables.clear();
372 string funName("SelectEVCalcVariables");
373 std::vector<std::string> variables = m_session->GetVariables();
374 for (size_t v = 0; v < m_nfields; ++v)
375 {
376 if (!m_session->DefinesFunction(funName, variables[v]))
377 {
378 unselectedVariables.insert(v);
379 }
380 }
381 if (unselectedVariables.size() == m_nfields)
382 {
383 unselectedVariables.clear();
384 }
385}
386
387/**
388 *
389 */
391{
392 std::vector<std::vector<LibUtilities::EquationSharedPtr>> selectedDomains;
393 std::set<int> unselectedVariables;
394 GetMaskInfo(selectedDomains, unselectedVariables);
395 if (selectedDomains.size() == 0 && unselectedVariables.size() == 0)
396 {
397 m_useMask = false;
398 return;
399 }
400 m_useMask = true;
401 MultiRegions::ExpListSharedPtr field = m_equ[0]->UpdateFields()[0];
402 int ncoef = field->GetNcoeffs();
403 int nphys = field->GetNpoints();
406 for (size_t i = 0; i < field->GetExpSize(); ++i)
407 {
408 LocalRegions::ExpansionSharedPtr exp = field->GetExp(i);
409 SpatialDomains::GeometrySharedPtr geom = exp->GetGeom();
410 int nv = geom->GetNumVerts();
411 NekDouble gc[3] = {0., 0., 0.};
412 NekDouble gct[3] = {0., 0., 0.};
413 for (size_t j = 0; j < nv; ++j)
414 {
415 SpatialDomains::PointGeomSharedPtr vertex = geom->GetVertex(j);
416 vertex->GetCoords(gct[0], gct[1], gct[2]);
417 gc[0] += gct[0] / NekDouble(nv);
418 gc[1] += gct[1] / NekDouble(nv);
419 gc[2] += gct[2] / NekDouble(nv);
420 }
421 int unmask = 1;
422 for (size_t m = 0; m < selectedDomains.size(); ++m)
423 {
424 unmask = 0;
425 for (size_t n = 0; n < selectedDomains[m].size(); ++n)
426 {
427 if (selectedDomains[m][n]->Evaluate(gc[0], gc[1], gc[2]) <= 0.)
428 {
429 unmask = 0;
430 break;
431 }
432 else
433 {
434 unmask = 1;
435 }
436 }
437 if (unmask == 1)
438 break;
439 }
440 for (int j = 0; j < m_nfields; ++j)
441 {
442 if (unmask == 0 || unselectedVariables.count(j))
443 {
445 exp->GetNcoeffs(), 0.,
446 &m_maskCoeffs[field->GetCoeff_Offset(i) + j * ncoef], 1);
447 Vmath::Fill(exp->GetTotPoints(), 0.,
448 &m_maskPhys[field->GetPhys_Offset(i) + j * nphys],
449 1);
450 }
451 }
452 }
453}
454
455} // namespace SolverUtils
456} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
void CopyFwdToAdj()
Copy the forward field to the adjoint system in transient growth calculations.
virtual 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:66
int m_infosteps
underlying operator is time stepping
Definition: DriverArnoldi.h:69
virtual void v_Execute(std::ostream &out=std::cout) override
Virtual function for solve implementation.
Array< OneD, NekDouble > m_maskCoeffs
Definition: DriverArnoldi.h:80
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:63
bool m_timeSteppingAlgorithm
Period of time stepping algorithm.
Definition: DriverArnoldi.h:67
int m_nits
Number of vectors to test.
Definition: DriverArnoldi.h:64
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:81
NekDouble m_evtol
Maxmum number of iterations.
Definition: DriverArnoldi.h:65
int m_nfields
interval to dump information if required.
Definition: DriverArnoldi.h:71
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:67
LibUtilities::SessionReaderSharedPtr m_session
Session reader object.
Definition: Driver.h:85
virtual SOLVER_UTILS_EXPORT void v_InitObject(std::ostream &out=std::cout)
Virtual function for initialisation implementation.
Definition: Driver.cpp:93
LibUtilities::CommSharedPtr m_comm
Communication object.
Definition: Driver.h:82
enum EvolutionOperatorType m_EvolutionOperator
Evolution Operator.
Definition: Driver.h:100
Array< OneD, EquationSystemSharedPtr > m_equ
Equation system to solve.
Definition: Driver.h:94
int m_nequ
number of equations
Definition: Driver.h:97
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kNekUnsetDouble
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:54
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303