Nektar++
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
42using namespace std;
43
44namespace Nektar
45{
46using namespace MultiRegions;
47
50 "SmoothedProfileMethod", SmoothedProfileMethod::create);
51
54 "SolverType", "SmoothedProfileMethod", eSmoothedProfileMethod);
55
56/**
57 * @brief Construct a new Smoothed Profile Method object
58 *
59 * @param pSession
60 * @param pGraph
61 */
65 : UnsteadySystem(pSession, pGraph),
66 VelocityCorrectionScheme(pSession, pGraph)
67{
68}
69
70/**
71 * @brief Destroy the Smoothed Profile Method object
72 *
73 */
75{
76}
77
79{
81
82 // Update implicit time-intregration class operators
85
86 // Number of dims as number of velocity vectors
87 size_t nvel = m_velocity.size();
88
89 // Initialization of correction pressure and shape function
90 switch (nvel)
91 {
92 case 1:
94 {
95 SetUpExpansions<ContField>(nvel);
96 }
98 {
99 SetUpExpansions<DisContField>(nvel);
100 }
101 break;
102
103 case 2:
105 {
106 SetUpExpansions<ContField>(nvel);
107 }
109 {
110 SetUpExpansions<DisContField>(nvel);
111 }
112 break;
113
114 case 3:
116 {
118 {
119 SetUpExpansions<ContField>(nvel);
120 }
122 {
123 SetUpExpansions<ContField3DHomogeneous1D>(nvel);
124 }
127 {
128 SetUpExpansions<ContField3DHomogeneous2D>(nvel);
129 }
130 }
132 {
134 {
135 SetUpExpansions<DisContField>(nvel);
136 }
138 {
139 SetUpExpansions<DisContField3DHomogeneous1D>(nvel);
140 }
143 {
144 SetUpExpansions<DisContField3DHomogeneous2D>(nvel);
145 }
146 }
147 break;
148 }
149
150 // Read 'm_phi' and its velocity
151 ASSERTL0(m_session->DefinesFunction("ShapeFunction"),
152 "ShapeFunction must be defined in the session file.")
153 ReadPhi();
154
155 // Allocate the vector 'm_up'
156 size_t physTot = m_phi->GetTotPoints();
157 m_velName.push_back("Up");
158 if (nvel > 1)
159 {
160 m_velName.push_back("Vp");
161 }
162 if (nvel == 3)
163 {
164 m_velName.push_back("Wp");
165 }
166
168 for (size_t i = 0; i < nvel; ++i)
169 {
170 m_up[i] = Array<OneD, NekDouble>(physTot, 0.0);
171 }
172
173 // Make sure that m_phi and m_up are defined
174 UpdatePhiUp(0.0);
175
176 // Get the time integration scheme.
178 if (m_session->DefinesTimeIntScheme())
179 {
180 timeInt = m_session->GetTimeIntScheme();
181 }
182 else
183 {
184 timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
185 timeInt.order = timeInt.method.back() - '0';
186
187 // Remove everything past the IMEX.
188 timeInt.method = timeInt.method.substr(0, 4);
189 }
190
191 // Select 'm_gamma0' depending on IMEX order
192 ASSERTL0(
193 boost::iequals(timeInt.method, "IMEX") && 1 <= timeInt.order &&
194 timeInt.order <= 4,
195 "The TimeIntegrationMethod scheme must be IMEX with order '1' to '4'.")
196
197 switch (timeInt.order)
198 {
199 case 1:
200 m_gamma0 = 1.0;
201 break;
202
203 case 2:
204 m_gamma0 = 3.0 / 2.0;
205 break;
206
207 case 3:
208 m_gamma0 = 11.0 / 6.0;
209 break;
210
211 case 4:
212 m_gamma0 = 25.0 / 12.0;
213 break;
214 }
215
216 // Check if the aeroforces filter is active, negative if inactive
217 m_forcesFilter = -1;
218 for (size_t i = 0; i < m_session->GetFilters().size(); ++i)
219 {
220 if (boost::iequals(m_session->GetFilters()[i].first, "AeroForcesSPM"))
221 {
222 m_forcesFilter = i;
223 break;
224 }
225 }
226}
227
228/**
229 * @brief Generates the summary of the current simulation
230 *
231 * @param s
232 */
234{
236 SolverUtils::AddSummaryItem(s, "IB formulation",
237 "Smoothed Profile Method (SPM)");
238}
239
240/**
241 * @brief Linear terms due to pressure and visosity are calculated here.
242 * After solving the velocity filed without taking into account the
243 * immersed boundaries, a new correction is applied through the force
244 * \f$f_s\f$:
245 *
246 * \f[ \mathbf{f_s} = \frac{\Phi^{n+1}(\mathbf{u_p}-\mathbf{u^*})}
247 * {\Delta t} \f]
248 *
249 * @param inarray
250 * @param outarray
251 * @param time
252 * @param a_iixDt
253 */
255 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
256 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
257 const NekDouble a_iixDt)
258{
260 time, a_iixDt);
261
262 size_t physTot = m_pressureP->GetNpoints();
263
264 /* SPM correction of velocity */
265 // Update 'm_phi' and 'm_up' if needed (evaluated at next time step)
266 UpdatePhiUp(time + a_iixDt);
267 // Update calculation of IB forcing 'm_fs'
268 UpdateForcing(outarray, a_iixDt);
269 // Estimate forces only if requested
270 if (m_forcesFilter >= 0)
271 {
272 static_pointer_cast<FilterAeroForcesSPM>(
274 ->CalculateForces(outarray, m_upPrev, m_phi, time, a_iixDt);
275 }
276 // Set BC conditions for pressure p_p
277 SetUpCorrectionPressure(outarray, m_F);
278 // Solve Poisson equation for pressure p_p
280 // Solve velocity in the next step with IB
281 SolveCorrectedVelocity(m_F, outarray, a_iixDt);
282
283 // Add pressures to get final value
284 Vmath::Vadd(physTot, m_pressure->GetPhys(), 1, m_pressureP->GetPhys(), 1,
285 m_pressure->UpdatePhys(), 1);
286 m_pressure->FwdTrans(m_pressure->GetPhys(), m_pressure->UpdateCoeffs());
287
288 // Add presure to outflow bc if using convective like BCs
289 m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
290}
291
292/**
293 * @brief Sets the forcing term of the equation for the correction pressure
294 * \f$p_p\f$:
295 *
296 * \f[ \nabla\cdot\mathbf{f_s} \f]
297 *
298 * @param fields
299 * @param Forcing
300 */
302 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &fields,
304{
305 size_t physTot = m_fs[0]->GetNpoints();
306 size_t nvel = m_velocity.size();
307
308 // Set boundary conditions
310
311 // Divergence of 'fs'
312 m_fields[m_velocity[0]]->PhysDeriv(eX, m_fs[0]->GetPhys(), Forcing[0]);
313
314 // Using 'Forcing[1]' as storage
315 for (size_t i = 1; i < nvel; ++i)
316 {
317 size_t ind = m_velocity[i];
318 m_fields[ind]->PhysDeriv(DirCartesianMap[i], m_fs[i]->GetPhys(),
319 Forcing[1]);
320 Vmath::Vadd(physTot, Forcing[1], 1, Forcing[0], 1, Forcing[0], 1);
321 }
322}
323
324/**
325 * @brief Solves the Poisson equation for the correction pressure
326 * \f$p_p\f$:
327 *
328 * \f[ \nabla^2 p_p = \nabla\cdot\mathbf{f_s} \f]
329 *
330 * @param Forcing
331 */
334{
336 // Factor 'lambda=0' in Helmholtz equation to get the Poisson form
338
339 // Solve the Poisson equation
340 m_pressureP->HelmSolve(Forcing, m_pressureP->UpdateCoeffs(), factors);
341
342 // Update node values from coefficients
343 m_pressureP->BwdTrans(m_pressureP->GetCoeffs(), m_pressureP->UpdatePhys());
344}
345
346/**
347 * @brief Corrects the velocity field so that the IBs are taken into
348 * account. Solves the explicit equation:
349 *
350 * \f[ \frac{\gamma_0(\mathbf{u_p}^{n+1} - \mathbf{u}^*)}{\Delta t} =
351 * \mathbf{f_s} - \nabla p_p \f]
352 *
353 * @param Forcing
354 * @param fields
355 * @param dt
356 */
360{
361 size_t physTot = m_phi->GetNpoints();
362
363 // Gradient of p_p
364 size_t nvel = m_velocity.size();
365 if (nvel == 2)
366 {
367 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1]);
368 }
369 else
370 {
371 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1],
372 Forcing[2]);
373 }
374
375 // Velocity correction
376 for (size_t i = 0; i < nvel; ++i)
377 {
378 // Adding -(1-m_phi)*grad(p_p) instead of -grad(p_p) reduces the
379 // flux through the walls, but the flow is not incompressible
380 if (m_session->DefinesSolverInfo("ForceBoundary") &&
381 boost::iequals(m_session->GetSolverInfo("ForceBoundary"), "True"))
382 {
383 Vmath::Vvtvm(physTot, m_phi->GetPhys(), 1, Forcing[i], 1,
384 Forcing[i], 1, Forcing[i], 1);
385 Vmath::Vadd(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
386 Forcing[i], 1);
387 }
388 else
389 {
390 Vmath::Vsub(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
391 Forcing[i], 1);
392 }
393 Blas::Daxpy(physTot, dt / m_gamma0, Forcing[i], 1, fields[i], 1);
394 }
395}
396
397/**
398 * @brief Updates the BCs for boundaries with Neumann or periodic BCs in
399 * the pressure:
400 *
401 * \f[ \frac{\partial p_p}{\partial\mathbf{n}} =
402 * \mathbf{f_s}\cdot\mathbf{n} \f]
403 */
405{
406 size_t nvel = m_velocity.size();
409
410 // Get the BC expansions
411 BndExp = m_pressureP->GetBndCondExpansions();
412 BndCond = m_pressureP->GetBndConditions();
413
414 // For each boundary...
415 for (size_t b = 0; b < BndExp.size(); ++b)
416 {
417 // Only for BCs based on the derivative
418 if (BndCond[b]->GetBoundaryConditionType() ==
420 BndCond[b]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
421 {
422 // Calculate f_s values
424 for (size_t i = 0; i < nvel; ++i)
425 {
426 f_s[i] = m_fs[0]->GetBndCondExpansions()[b]->GetPhys();
427 }
428
429 // BC is f_s * n
430 BndExp[b]->NormVectorIProductWRTBase(f_s, BndExp[b]->UpdatePhys());
431 }
432 }
433}
434
435/**
436 * @brief Calculates the values of the shape function
437 *
438 * @param t
439 */
441{
442 // Initialise 'm_up' and 'm_phi' during first step
443 if (t <= 0.0)
444 {
445 if (!m_filePhi)
446 {
447 // Update 'm_phi' only if it was provided as a function
448 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
449 }
450
451 // Initialize both variables for the first step
452 m_phiEvaluator->Evaluate(m_velName, m_up, t);
453
454 // Initialise 'm_upPrev' in all cases
455 m_upPrev = m_up;
456 }
457 // If timedependent 'm_phi'
458 // Phi functions from files are not timedependent
459 else if (m_timeDependentPhi)
460 {
461 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
462
463 // And if velocities are timedependent as well
465 {
466 // Store previous value of u_p during simulation
467 m_upPrev = m_up;
468 m_phiEvaluator->Evaluate(m_velName, m_up, t);
469 }
470 }
471}
472
473/**
474 * @brief For a body with a velocity \f$\mathbf{u_p}\f$, the force
475 * \f$\mathbf{f_s}\f$ applied to the fluid ensures that the IBC are met:
476 *
477 * \f[ \mathbf{f_s} = \frac{\Phi^{n+1}\left(\mathbf{u_p}^{n+1} -
478 * \mathbf{u^*}\right)}{\Delta t} \f]
479 *
480 * @param fields
481 * @param dt
482 * @param f_s
483 */
485 const Array<OneD, const Array<OneD, NekDouble>> &fields, NekDouble dt)
486{
487 size_t nvel = m_velocity.size();
488 size_t nq = m_phi->GetNpoints();
489
490 for (size_t i = 0; i < nvel; ++i)
491 {
492 // In homogeneous cases, switch out of wave space
493 Array<OneD, NekDouble> tmpField(nq);
494 size_t ind = m_velocity[i];
495
497 m_fields[ind]->GetWaveSpace())
498 {
499 m_fields[ind]->HomogeneousBwdTrans(nq, fields[i], tmpField);
500 m_fs[i]->HomogeneousBwdTrans(nq, m_fs[i]->GetPhys(),
501 m_fs[i]->UpdatePhys());
502 }
503 else
504 {
505 tmpField = fields[i];
506 }
507
508 Vmath::Vsub(nq, m_up[i], 1, tmpField, 1, m_fs[i]->UpdatePhys(), 1);
509 Vmath::Vmul(nq, m_phi->GetPhys(), 1, m_fs[i]->GetPhys(), 1,
510 m_fs[i]->UpdatePhys(), 1);
511 Vmath::Smul(nq, m_gamma0 / dt, m_fs[i]->GetPhys(), 1,
512 m_fs[i]->UpdatePhys(), 1);
513
514 // And go back to wave space if the 'if' was executed
516 m_fields[ind]->GetWaveSpace())
517 {
518 m_fs[i]->HomogeneousFwdTrans(nq, m_fs[i]->GetPhys(),
519 m_fs[i]->UpdatePhys());
520 }
521 }
522}
523
524/**
525 * @brief True if the function is timedependent, false otherwise
526 *
527 * @param name
528 * @param type
529 * @param attribute
530 * @return string
531 */
533 string elemName)
534{
535 // Get the handler of the function
536 TiXmlElement *function = GetFunctionHdl(funcName);
537
538 // Go to the first element
539 TiXmlElement *functionDef = function->FirstChildElement();
540 ASSERTL0(functionDef, "At least one element must be defined in " + funcName)
541
542 // And search the element with name 'elemName'
543 string varName = functionDef->Attribute("VAR");
544 while (functionDef && !boost::iequals(varName, elemName))
545 {
546 functionDef = functionDef->NextSiblingElement();
547 varName = functionDef->Attribute("VAR");
548 }
549
550 ASSERTL0(functionDef,
551 "Variable " + elemName + " must be defined in " + funcName + ".");
552
553 // And return the value of USERDEFINEDTYPE
554 string attr;
555 int err = functionDef->QueryStringAttribute("USERDEFINEDTYPE", &attr);
556 bool output = boost::iequals(attr, "TimeDependent");
557
558 ASSERTL0((err == TIXML_NO_ATTRIBUTE) || (err == TIXML_SUCCESS && output),
559 "USERDEFINEDTYPE in " + elemName +
560 " must be TimeDependent if defined");
561
562 return output;
563}
564
565/**
566 * @brief Returns a handle to the requested function. Returns NULL if it
567 * does not exist
568 *
569 * @param functionName
570 * @return TiXmlElement*
571 */
572TiXmlElement *SmoothedProfileMethod::GetFunctionHdl(string functionName)
573{
574 // Get the handler of first function block
575 TiXmlElement *conds = m_session->GetElement("Nektar/Conditions");
576 TiXmlElement *function = conds->FirstChildElement("FUNCTION");
577
578 // Loop over functions until the block 'name' is found
579 string functionType = function->Attribute("NAME");
580 while (function && !boost::iequals(functionType, functionName))
581 {
582 function = function->NextSiblingElement("FUNCTION");
583 functionType = function->Attribute("NAME");
584 }
585
586 return function;
587}
588
590{
591 // Function evaluator for Phi and Up
592 m_phiEvaluator = GetFunction("ShapeFunction");
593
594 TiXmlElement *function = GetFunctionHdl("ShapeFunction");
595 TiXmlElement *child = function->FirstChildElement();
596 m_filePhi = false;
597
598 // If defined by using a file
599 if (boost::iequals(child->ValueStr(), "F"))
600 {
601 // Get name of STL file
602 string fileName;
603 int status = child->QueryStringAttribute("FILE", &fileName);
604 ASSERTL0(status == TIXML_SUCCESS,
605 "An FLD file with the values "
606 "of the phi function has to be supplied.")
607 ASSERTL0(boost::iequals(fileName.substr(fileName.length() - 4), ".fld"),
608 "A valid FLD file must be supplied in the "
609 "'ShapeFunction' field.")
610
611 // Get phi values from XML file (after "FieldConvert" the STL file)
612 // First, load the data
613 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
614 std::vector<std::vector<NekDouble>> fieldData;
615 LibUtilities::FieldMetaDataMap fieldMetaData;
618 phiFile->Import(fileName, fieldDef, fieldData, fieldMetaData);
619
620 // Only Phi field should be defined in the file
621 ASSERTL0(fieldData.size() == 1, "Only one field (phi) must be "
622 "defined in the FLD file.")
623
624 // Extract Phi field to output
625 string tmp("phi");
626 m_phi->ExtractDataToCoeffs(fieldDef[0], fieldData[0], tmp,
627 m_phi->UpdateCoeffs());
628 m_phi->BwdTrans(m_phi->GetCoeffs(), m_phi->UpdatePhys());
629 m_filePhi = true;
630 m_timeDependentPhi = false;
631 m_timeDependentUp = false;
632 }
633 else
634 {
635 // Check if Phi is timedependent
636 m_timeDependentPhi = GetVarTimeDependence("ShapeFunction", "Phi");
637
638 // If so, check if its velocity changes as well
639 m_timeDependentUp = GetVarTimeDependence("ShapeFunction", "Up");
640 switch (m_velocity.size())
641 {
642 case 2:
644 GetVarTimeDependence("ShapeFunction", "Vp");
645 break;
646 case 3:
648 GetVarTimeDependence("ShapeFunction", "Vp");
650 GetVarTimeDependence("ShapeFunction", "Wp");
651 break;
652 }
653 }
654}
655
656} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:224
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.
~SmoothedProfileMethod() override
Destroy the Smoothed Profile Method object.
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 HomogeneousType m_HomogeneousType
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:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
@ eSmoothedProfileMethod
double NekDouble
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
STL namespace.