Nektar++
Loading...
Searching...
No Matches
SmoothedProfileMethod.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: SmoothedProfileMethod.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// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Smoothed Profile Method for the Incompressible
33// Navier Stokes equations
34///////////////////////////////////////////////////////////////////////////////
35
41
42namespace Nektar
43{
44
45using namespace MultiRegions;
46
49 "SmoothedProfileMethod", SmoothedProfileMethod::create);
50
53 "SolverType", "SmoothedProfileMethod", eSmoothedProfileMethod);
54
55/**
56 * @brief Construct a new Smoothed Profile Method object
57 *
58 * @param pSession
59 * @param pGraph
60 */
68
70{
72
73 // Update implicit time-intregration class operators
76
77 // Number of dims as number of velocity vectors
78 size_t nvel = m_velocity.size();
79
80 // Initialization of correction pressure and shape function
81 switch (nvel)
82 {
83 case 1:
85 {
86 SetUpExpansions<ContField>(nvel);
87 }
89 {
90 SetUpExpansions<DisContField>(nvel);
91 }
92 break;
93
94 case 2:
96 {
97 SetUpExpansions<ContField>(nvel);
98 }
100 {
101 SetUpExpansions<DisContField>(nvel);
102 }
103 break;
104
105 case 3:
107 {
109 {
110 SetUpExpansions<ContField>(nvel);
111 }
113 {
114 SetUpExpansions<ContField3DHomogeneous1D>(nvel);
115 }
118 {
119 SetUpExpansions<ContField3DHomogeneous2D>(nvel);
120 }
121 }
123 {
125 {
126 SetUpExpansions<DisContField>(nvel);
127 }
129 {
130 SetUpExpansions<DisContField3DHomogeneous1D>(nvel);
131 }
134 {
135 SetUpExpansions<DisContField3DHomogeneous2D>(nvel);
136 }
137 }
138 break;
139 }
140
141 // Read 'm_phi' and its velocity
142 ASSERTL0(m_session->DefinesFunction("ShapeFunction"),
143 "ShapeFunction must be defined in the session file.")
144 ReadPhi();
145
146 // Allocate the vector 'm_up'
147 size_t physTot = m_phi->GetTotPoints();
148 m_velName.push_back("Up");
149 if (nvel > 1)
150 {
151 m_velName.push_back("Vp");
152 }
153 if (nvel == 3)
154 {
155 m_velName.push_back("Wp");
156 }
157
159 for (size_t i = 0; i < nvel; ++i)
160 {
161 m_up[i] = Array<OneD, NekDouble>(physTot, 0.0);
162 }
163
164 // Make sure that m_phi and m_up are defined
165 UpdatePhiUp(0.0);
166
167 // Get the time integration scheme.
168 LibUtilities::TimeIntScheme timeInt = m_session->GetTimeIntScheme();
169
170 // Select 'm_gamma0' depending on IMEX order
171 ASSERTL0(
172 boost::iequals(timeInt.method, "IMEX") && 1 <= timeInt.order &&
173 timeInt.order <= 4,
174 "The TimeIntegrationMethod scheme must be IMEX with order '1' to '4'.")
175
176 switch (timeInt.order)
177 {
178 case 1:
179 m_gamma0 = 1.0;
180 break;
181
182 case 2:
183 m_gamma0 = 3.0 / 2.0;
184 break;
185
186 case 3:
187 m_gamma0 = 11.0 / 6.0;
188 break;
189
190 case 4:
191 m_gamma0 = 25.0 / 12.0;
192 break;
193 }
194
195 // Check if the aeroforces filter is active, negative if inactive
196 m_forcesFilter = -1;
197 for (size_t i = 0; i < m_session->GetFilters().size(); ++i)
198 {
199 if (boost::iequals(m_session->GetFilters()[i].name, "AeroForcesSPM"))
200 {
201 m_forcesFilter = i;
202 break;
203 }
204 }
205}
206
207/**
208 * @brief Generates the summary of the current simulation
209 *
210 * @param s
211 */
213{
215 SolverUtils::AddSummaryItem(s, "IB formulation",
216 "Smoothed Profile Method (SPM)");
217}
218
219/**
220 * @brief Linear terms due to pressure and visosity are calculated here.
221 * After solving the velocity filed without taking into account the
222 * immersed boundaries, a new correction is applied through the force
223 * \f$f_s\f$:
224 *
225 * \f[ \mathbf{f_s} = \frac{\Phi^{n+1}(\mathbf{u_p}-\mathbf{u^*})}
226 * {\Delta t} \f]
227 *
228 * @param inarray
229 * @param outarray
230 * @param time
231 * @param a_iixDt
232 */
234 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
235 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
236 const NekDouble a_iixDt)
237{
239 time, a_iixDt);
240
241 size_t physTot = m_pressureP->GetNpoints();
242
243 /* SPM correction of velocity */
244 // Update 'm_phi' and 'm_up' if needed (evaluated at next time step)
245 UpdatePhiUp(time + a_iixDt);
246 // Update calculation of IB forcing 'm_fs'
247 UpdateForcing(outarray, a_iixDt);
248 // Estimate forces only if requested
249 if (m_forcesFilter >= 0)
250 {
251 std::static_pointer_cast<FilterAeroForcesSPM>(
253 ->CalculateForces(outarray, m_upPrev, m_phi, time, a_iixDt);
254 }
255 // Set BC conditions for pressure p_p
256 SetUpCorrectionPressure(outarray, m_F);
257 // Solve Poisson equation for pressure p_p
259 // Solve velocity in the next step with IB
260 SolveCorrectedVelocity(m_F, outarray, a_iixDt);
261
262 // Add pressures to get final value
263 Vmath::Vadd(physTot, m_pressure->GetPhys(), 1, m_pressureP->GetPhys(), 1,
264 m_pressure->UpdatePhys(), 1);
265 m_pressure->FwdTrans(m_pressure->GetPhys(), m_pressure->UpdateCoeffs());
266
267 // Add presure to outflow bc if using convective like BCs
268 m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
269}
270
271/**
272 * @brief Sets the forcing term of the equation for the correction pressure
273 * \f$p_p\f$:
274 *
275 * \f[ \nabla\cdot\mathbf{f_s} \f]
276 *
277 * @param fields
278 * @param Forcing
279 */
281 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &fields,
283{
284 size_t physTot = m_fs[0]->GetNpoints();
285 size_t nvel = m_velocity.size();
286
287 // Set boundary conditions
289
290 // Divergence of 'fs'
291 m_fields[m_velocity[0]]->PhysDeriv(eX, m_fs[0]->GetPhys(), Forcing[0]);
292
293 // Using 'Forcing[1]' as storage
294 for (size_t i = 1; i < nvel; ++i)
295 {
296 size_t ind = m_velocity[i];
297 m_fields[ind]->PhysDeriv(DirCartesianMap[i], m_fs[i]->GetPhys(),
298 Forcing[1]);
299 Vmath::Vadd(physTot, Forcing[1], 1, Forcing[0], 1, Forcing[0], 1);
300 }
301}
302
303/**
304 * @brief Solves the Poisson equation for the correction pressure
305 * \f$p_p\f$:
306 *
307 * \f[ \nabla^2 p_p = \nabla\cdot\mathbf{f_s} \f]
308 *
309 * @param Forcing
310 */
313{
315 // Factor 'lambda=0' in Helmholtz equation to get the Poisson form
316 factors[StdRegions::eFactorLambda] = 0.0;
317
318 // Solve the Poisson equation
319 m_pressureP->HelmSolve(Forcing, m_pressureP->UpdateCoeffs(), factors);
320
321 // Update node values from coefficients
322 m_pressureP->BwdTrans(m_pressureP->GetCoeffs(), m_pressureP->UpdatePhys());
323}
324
325/**
326 * @brief Corrects the velocity field so that the IBs are taken into
327 * account. Solves the explicit equation:
328 *
329 * \f[ \frac{\gamma_0(\mathbf{u_p}^{n+1} - \mathbf{u}^*)}{\Delta t} =
330 * \mathbf{f_s} - \nabla p_p \f]
331 *
332 * @param Forcing
333 * @param fields
334 * @param dt
335 */
339{
340 size_t physTot = m_phi->GetNpoints();
341
342 // Gradient of p_p
343 size_t nvel = m_velocity.size();
344 if (nvel == 2)
345 {
346 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1]);
347 }
348 else
349 {
350 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1],
351 Forcing[2]);
352 }
353
354 // Velocity correction
355 for (size_t i = 0; i < nvel; ++i)
356 {
357 // Adding -(1-m_phi)*grad(p_p) instead of -grad(p_p) reduces the
358 // flux through the walls, but the flow is not incompressible
359 if (m_session->DefinesSolverInfo("ForceBoundary") &&
360 boost::iequals(m_session->GetSolverInfo("ForceBoundary"), "True"))
361 {
362 Vmath::Vvtvm(physTot, m_phi->GetPhys(), 1, Forcing[i], 1,
363 Forcing[i], 1, Forcing[i], 1);
364 Vmath::Vadd(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
365 Forcing[i], 1);
366 }
367 else
368 {
369 Vmath::Vsub(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
370 Forcing[i], 1);
371 }
372 Blas::Daxpy(physTot, dt / m_gamma0, Forcing[i], 1, fields[i], 1);
373 }
374}
375
376/**
377 * @brief Updates the BCs for boundaries with Neumann or periodic BCs in
378 * the pressure:
379 *
380 * \f[ \frac{\partial p_p}{\partial\mathbf{n}} =
381 * \mathbf{f_s}\cdot\mathbf{n} \f]
382 */
384{
385 size_t nvel = m_velocity.size();
388
389 // Get the BC expansions
390 BndExp = m_pressureP->GetBndCondExpansions();
391 BndCond = m_pressureP->GetBndConditions();
392
393 // For each boundary...
394 for (size_t b = 0; b < BndExp.size(); ++b)
395 {
396 // Only for BCs based on the derivative
397 if (BndCond[b]->GetBoundaryConditionType() ==
399 BndCond[b]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
400 {
401 // Calculate f_s values
403 for (size_t i = 0; i < nvel; ++i)
404 {
405 f_s[i] = m_fs[0]->GetBndCondExpansions()[b]->GetPhys();
406 }
407
408 // BC is f_s * n
409 BndExp[b]->NormVectorIProductWRTBase(f_s, BndExp[b]->UpdatePhys());
410 }
411 }
412}
413
414/**
415 * @brief Calculates the values of the shape function
416 *
417 * @param t
418 */
420{
421 // Initialise 'm_up' and 'm_phi' during first step
422 if (t <= 0.0)
423 {
424 if (!m_filePhi)
425 {
426 // Update 'm_phi' only if it was provided as a function
427 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
428 }
429
430 // Initialize both variables for the first step
431 m_phiEvaluator->Evaluate(m_velName, m_up, t);
432
433 // Initialise 'm_upPrev' in all cases
434 m_upPrev = m_up;
435 }
436 // If timedependent 'm_phi'
437 // Phi functions from files are not timedependent
438 else if (m_timeDependentPhi)
439 {
440 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
441
442 // And if velocities are timedependent as well
444 {
445 // Store previous value of u_p during simulation
446 m_upPrev = m_up;
447 m_phiEvaluator->Evaluate(m_velName, m_up, t);
448 }
449 }
450}
451
452/**
453 * @brief For a body with a velocity \f$\mathbf{u_p}\f$, the force
454 * \f$\mathbf{f_s}\f$ applied to the fluid ensures that the IBC are met:
455 *
456 * \f[ \mathbf{f_s} = \frac{\Phi^{n+1}\left(\mathbf{u_p}^{n+1} -
457 * \mathbf{u^*}\right)}{\Delta t} \f]
458 *
459 * @param fields
460 * @param dt
461 * @param f_s
462 */
464 const Array<OneD, const Array<OneD, NekDouble>> &fields, NekDouble dt)
465{
466 size_t nvel = m_velocity.size();
467 size_t nq = m_phi->GetNpoints();
468
469 for (size_t i = 0; i < nvel; ++i)
470 {
471 // In homogeneous cases, switch out of wave space
472 Array<OneD, NekDouble> tmpField(nq);
473 size_t ind = m_velocity[i];
474
476 m_fields[ind]->GetWaveSpace())
477 {
478 m_fields[ind]->HomogeneousBwdTrans(nq, fields[i], tmpField);
479 m_fs[i]->HomogeneousBwdTrans(nq, m_fs[i]->GetPhys(),
480 m_fs[i]->UpdatePhys());
481 }
482 else
483 {
484 tmpField = fields[i];
485 }
486
487 Vmath::Vsub(nq, m_up[i], 1, tmpField, 1, m_fs[i]->UpdatePhys(), 1);
488 Vmath::Vmul(nq, m_phi->GetPhys(), 1, m_fs[i]->GetPhys(), 1,
489 m_fs[i]->UpdatePhys(), 1);
490 Vmath::Smul(nq, m_gamma0 / dt, m_fs[i]->GetPhys(), 1,
491 m_fs[i]->UpdatePhys(), 1);
492
493 // And go back to wave space if the 'if' was executed
495 m_fields[ind]->GetWaveSpace())
496 {
497 m_fs[i]->HomogeneousFwdTrans(nq, m_fs[i]->GetPhys(),
498 m_fs[i]->UpdatePhys());
499 }
500 }
501}
502
503/**
504 * @brief True if the function is timedependent, false otherwise
505 *
506 * @param name
507 * @param type
508 * @param attribute
509 * @return string
510 */
512 std::string elemName)
513{
514 // Get the handler of the function
515 TiXmlElement *function = GetFunctionHdl(funcName);
516
517 // Go to the first element
518 TiXmlElement *functionDef = function->FirstChildElement();
519 ASSERTL0(functionDef, "At least one element must be defined in " + funcName)
520
521 // And search the element with name 'elemName'
522 std::string varName = functionDef->Attribute("VAR");
523 while (functionDef && !boost::iequals(varName, elemName))
524 {
525 functionDef = functionDef->NextSiblingElement();
526 varName = functionDef->Attribute("VAR");
527 }
528
529 ASSERTL0(functionDef,
530 "Variable " + elemName + " must be defined in " + funcName + ".");
531
532 // And return the value of USERDEFINEDTYPE
533 std::string attr;
534 int err = functionDef->QueryStringAttribute("USERDEFINEDTYPE", &attr);
535 bool output = boost::iequals(attr, "TimeDependent");
536
537 ASSERTL0((err == TIXML_NO_ATTRIBUTE) || (err == TIXML_SUCCESS && output),
538 "USERDEFINEDTYPE in " + elemName +
539 " must be TimeDependent if defined");
540
541 return output;
542}
543
544/**
545 * @brief Returns a handle to the requested function. Returns NULL if it
546 * does not exist
547 *
548 * @param functionName
549 * @return TiXmlElement*
550 */
551TiXmlElement *SmoothedProfileMethod::GetFunctionHdl(std::string functionName)
552{
553 // Get the handler of first function block
554 TiXmlElement *conds = m_session->GetElement("Nektar/Conditions");
555 TiXmlElement *function = conds->FirstChildElement("FUNCTION");
556
557 // Loop over functions until the block 'name' is found
558 std::string functionType = function->Attribute("NAME");
559 while (function && !boost::iequals(functionType, functionName))
560 {
561 function = function->NextSiblingElement("FUNCTION");
562 functionType = function->Attribute("NAME");
563 }
564
565 return function;
566}
567
569{
570 // Function evaluator for Phi and Up
571 m_phiEvaluator = GetFunction("ShapeFunction");
572
573 TiXmlElement *function = GetFunctionHdl("ShapeFunction");
574 TiXmlElement *child = function->FirstChildElement();
575 m_filePhi = false;
576
577 // If defined by using a file
578 if (boost::iequals(child->ValueStr(), "F"))
579 {
580 // Get name of STL file
581 std::string fileName;
582 int status = child->QueryStringAttribute("FILE", &fileName);
583 ASSERTL0(status == TIXML_SUCCESS,
584 "An FLD file with the values "
585 "of the phi function has to be supplied.")
586 ASSERTL0(boost::iequals(fileName.substr(fileName.length() - 4), ".fld"),
587 "A valid FLD file must be supplied in the "
588 "'ShapeFunction' field.")
589
590 // Get phi values from XML file (after "FieldConvert" the STL file)
591 // First, load the data
592 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
593 std::vector<std::vector<NekDouble>> fieldData;
594 LibUtilities::FieldMetaDataMap fieldMetaData;
597 phiFile->Import(fileName, fieldDef, fieldData, fieldMetaData);
598
599 // Only Phi field should be defined in the file
600 ASSERTL0(fieldData.size() == 1, "Only one field (phi) must be "
601 "defined in the FLD file.")
602
603 // Extract Phi field to output
604 std::string tmp("phi");
605 m_phi->ExtractDataToCoeffs(fieldDef[0], fieldData[0], tmp,
606 m_phi->UpdateCoeffs());
607 m_phi->BwdTrans(m_phi->GetCoeffs(), m_phi->UpdatePhys());
608 m_filePhi = true;
609 m_timeDependentPhi = false;
610 m_timeDependentUp = false;
611 }
612 else
613 {
614 // Check if Phi is timedependent
615 m_timeDependentPhi = GetVarTimeDependence("ShapeFunction", "Phi");
616
617 // If so, check if its velocity changes as well
618 m_timeDependentUp = GetVarTimeDependence("ShapeFunction", "Up");
619 switch (m_velocity.size())
620 {
621 case 2:
623 GetVarTimeDependence("ShapeFunction", "Vp");
624 break;
625 case 3:
627 GetVarTimeDependence("ShapeFunction", "Vp");
629 GetVarTimeDependence("ShapeFunction", "Wp");
630 break;
631 }
632 }
633}
634
635} // namespace Nektar
#define ASSERTL0(condition, msg)
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition FieldIO.cpp:223
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void UpdateForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble dt)
For a body with a velocity , the force applied to the fluid ensures that the IBC are met:
bool m_filePhi
Flag indicating that phi was defined in a file.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Generates the summary of the current simulation.
void SetCorrectionPressureBCs()
Updates the BCs for boundaries with Neumann or periodic BCs in the pressure:
void v_SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble a_iixDt) override
Linear terms due to pressure and visosity are calculated here. After solving the velocity filed witho...
Array< OneD, Array< OneD, NekDouble > > m_up
Velocity of the immersed body(ies)
std::vector< std::string > m_velName
Vector storing the names of the components of \u_p.
NekDouble m_gamma0
Stiffly-stable scheme coefficient.
int m_forcesFilter
Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'.
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
bool m_timeDependentPhi
Flag that is true when phi depends on time.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
SmoothedProfileMethod(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Construct a new Smoothed Profile Method object.
void UpdatePhiUp(NekDouble time)
Calculates the values of the shape function.
SolverUtils::SessionFunctionSharedPtr m_phiEvaluator
Function that evaluates the values of \Phi.
void SolveCorrectedVelocity(Array< OneD, Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble dt)
Corrects the velocity field so that the IBs are taken into account. Solves the explicit equation:
TiXmlElement * GetFunctionHdl(std::string functionName)
Returns a handle to the requested function. Returns NULL if it does not exist.
static std::string className
Name of class.
bool m_timeDependentUp
Flag signaling if depends on time.
MultiRegions::ExpListSharedPtr m_pressureP
Correction pressure field for SPM.
bool GetVarTimeDependence(std::string funcName, std::string attrName)
True if the function is timedependent, false otherwise.
Array< OneD, Array< OneD, NekDouble > > m_upPrev
Array< OneD, MultiRegions::ExpListSharedPtr > m_fs
Forcing function 'f_s'.
MultiRegions::ExpListSharedPtr m_phi
Shape function 'phi' as expansion list.
void SolveCorrectionPressure(const Array< OneD, NekDouble > &Forcing)
Solves the Poisson equation for the correction pressure :
void SetUpCorrectionPressure(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing)
Sets the forcing term of the equation for the correction pressure :
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
Defines a forcing term to be explicitly applied.
Definition Forcing.h:71
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, Array< OneD, NekDouble > > m_F
virtual void v_SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition Blas.hpp:135
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition FieldIO.h:50
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition ExpList.h:87
std::vector< std::pair< std::string, std::string > > SummaryList
Definition Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition Misc.cpp:47
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition MeshGraph.h:217
std::map< ConstFactorType, NekDouble > ConstFactorMap
@ eSmoothedProfileMethod
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 Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition Vmath.hpp:180
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
Definition Vmath.hpp:381
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 Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition Vmath.hpp:220