LagrangianDdtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LagrangianDdtScheme.H"
27 #include "LagrangianSubFields.H"
28 
29 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
30 
31 template<class Type>
34 (
35  const LagrangianMesh& mesh,
36  Istream& is
37 )
38 {
39  if (is.eof())
40  {
42  << "Ddt scheme not specified" << endl << endl
43  << "Valid ddt schemes are :" << endl
44  << IstreamConstructorTablePtr_->sortedToc()
45  << exit(FatalIOError);
46  }
47 
48  const word schemeName(is);
49 
50  typename IstreamConstructorTable::iterator cstrIter =
51  IstreamConstructorTablePtr_->find(schemeName);
52 
53  if (cstrIter == IstreamConstructorTablePtr_->end())
54  {
56  << "Unknown ddt scheme " << schemeName << nl << nl
57  << "Valid ddt schemes are :" << endl
58  << IstreamConstructorTablePtr_->sortedToc()
59  << exit(FatalIOError);
60  }
61 
62  return cstrIter()(mesh, is);
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
68 template<class Type>
70 {}
71 
72 
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
74 
75 template<class Type>
78 (
79  const LagrangianSubScalarField& deltaT,
81 )
82 {
84 
85  tEqn.ref().deltaTSp += 1;
86  tEqn.ref().deltaTSu -= psi.oldTime();
87 
88  return tEqn;
89 }
90 
91 
92 template<class Type>
95 (
96  const LagrangianSubScalarField& deltaT,
99 )
100 {
101  tmp<LagrangianEqn<Type>> tEqn(new LagrangianEqn<Type>(deltaT, psi));
102 
103  tEqn.ref().deltaTSp += m;
104  tEqn.ref().deltaTSu -= m.oldTime()*psi.oldTime();
105 
106  return tEqn;
107 }
108 
109 
110 template<class Type>
113 (
114  const LagrangianSubScalarField& deltaT,
116 )
117 {
118  tmp<LagrangianEqn<Type>> tEqn(new LagrangianEqn<Type>(deltaT, psi));
119 
120  tEqn.ref().deltaTSu += psi - psi.oldTime();
121 
122  return tEqn;
123 }
124 
125 
126 template<class Type>
129 (
130  const LagrangianSubScalarField& deltaT,
133 )
134 {
135  tmp<LagrangianEqn<Type>> tEqn(new LagrangianEqn<Type>(deltaT, psi));
136 
137  tEqn.ref().deltaTSu += m*psi - m.oldTime()*psi.oldTime();
138 
139  return tEqn;
140 }
141 
142 
143 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
Class containing Lagrangian geometry and topology.
virtual ~ddtScheme()
Destructor.
static tmp< ddtScheme< Type > > New(const LagrangianMesh &mesh, Istream &is)
Return a pointer to a new ddtScheme.
static tmp< LagrangianEqn< Type > > Lagrangianmddt0(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &psi)
Return the explicit/forward time-derivative matrix.
static tmp< LagrangianEqn< Type > > Lagrangianmddt(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &psi)
Return the time-derivative matrix.
const Field0Type & oldTime() const
Return the old-time field.
Definition: OldTimeField.C:322
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const volScalarField & psi
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
IOerror FatalIOError
static const char nl
Definition: Ostream.H:267