timeVaryingMassSorptionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "EulerDdtScheme.H"
33 #include "CrankNicolsonDdtScheme.H"
34 #include "backwardDdtScheme.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 const Foam::Enum
39 <
41 >
42 Foam::timeVaryingMassSorptionFvPatchScalarField::ddtSchemeTypeNames_
43 ({
44  {
45  ddtSchemeType::tsEuler,
46  fv::EulerDdtScheme<scalar>::typeName_()
47  },
48  {
49  ddtSchemeType::tsCrankNicolson,
50  fv::CrankNicolsonDdtScheme<scalar>::typeName_()
51  },
52  {
53  ddtSchemeType::tsBackward,
54  fv::backwardDdtScheme<scalar>::typeName_()
55  },
56 });
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
63 (
64  const fvPatch& p,
66 )
67 :
68  fixedValueFvPatchScalarField(p, iF),
69  kabs_(scalar(1)),
70  max_(scalar(1)),
71  kdes_(scalar(1))
72 {}
73 
74 
77 (
78  const fvPatch& p,
80  const dictionary& dict
81 )
82 :
83  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
84  kabs_(dict.getCheck<scalar>("kabs", scalarMinMax::ge(0))),
85  max_(dict.getCheck<scalar>("max", scalarMinMax::ge(0))),
86  kdes_(dict.getCheckOrDefault<scalar>("kdes", 0, scalarMinMax::ge(0)))
87 {
88  if (!this->readValueEntry(dict))
89  {
91  }
92 }
93 
94 
97 (
98  const timeVaryingMassSorptionFvPatchScalarField& ptf,
99  const fvPatch& p,
100  const DimensionedField<scalar, volMesh>& iF,
101  const fvPatchFieldMapper& mapper
102 )
103 :
104  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
105  kabs_(ptf.kabs_),
106  max_(ptf.max_),
107  kdes_(ptf.kdes_)
108 {}
109 
110 
113 (
115 )
116 :
117  fixedValueFvPatchScalarField(ptf),
118  kabs_(ptf.kabs_),
119  max_(ptf.max_),
120  kdes_(ptf.kdes_)
121 {}
122 
123 
126 (
129 )
130 :
131  fixedValueFvPatchScalarField(ptf, iF),
132  kabs_(ptf.kabs_),
133  max_(ptf.max_),
134  kdes_(ptf.kdes_)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 (
142  const fvPatchFieldMapper& m
143 )
144 {
145  fixedValueFvPatchScalarField::autoMap(m);
146 }
147 
148 
150 (
151  const fvPatchScalarField& ptf,
152  const labelList& addr
153 )
154 {
155  fixedValueFvPatchScalarField::rmap(ptf, addr);
156 }
157 
158 
161 {
162  auto tsource = tmp<scalarField>::New(patch().size(), Zero);
163  auto& source = tsource.ref();
164 
165  const scalarField cp(*this);
166  const scalarField w(max(1 - cp/max_, scalar(0)));
167 
168  source = -kabs_*w*max(patchInternalField() - cp, scalar(0));
170  source += kdes_*max(cp - patchInternalField(), scalar(0));
171 
172  return tsource;
173 }
174 
175 
177 {
178  if (updated())
179  {
180  return;
181  }
182 
183  const label patchi = patch().index();
184 
185  const scalar dt = db().time().deltaTValue();
186 
187  const auto& fld =
188  db().lookupObject<volScalarField>(this->internalField().name());
189  const volScalarField& fld0 = fld.oldTime();
190 
191  // Lookup d/dt scheme from database
192  const word ddtSchemeName(fld.mesh().ddtScheme(fld.name()));
193  const ddtSchemeType& ddtScheme = ddtSchemeTypeNames_[ddtSchemeName];
194 
195  const scalarField cp(*this);
196  const scalarField w(max(1 - cp/max_, scalar(0)));
197 
198  tmp<scalarField> dfldp =
199  kabs_
200  *w
201  *max(patchInternalField() - cp, scalar(0))
202  *dt;
203 
204  dfldp.ref() -=
205  kdes_
206  *max(cp - patchInternalField(), scalar(0))
207  *dt;
208 
209  switch (ddtScheme)
210  {
211  case tsEuler:
212  case tsCrankNicolson:
213  {
214  operator==(fld0.boundaryField()[patchi] + dfldp);
215 
216  break;
217  }
218  case tsBackward:
219  {
220  const scalar dt0 = db().time().deltaT0Value();
221 
222  const scalar c = scalar(1) + dt/(dt + dt0);
223  const scalar c00 = dt*dt/(dt0*(dt + dt0));
224  const scalar c0 = c + c00;
225 
226  operator==
227  (
228  (
229  c0*fld0.boundaryField()[patchi]
230  - c00*fld0.oldTime().boundaryField()[patchi]
231  + dfldp
232  )/c
233  );
234 
235  break;
236  }
237  default:
238  {
240  << ddtSchemeName << nl
241  << " on patch " << this->patch().name()
242  << " of field " << this->internalField().name()
243  << " in file " << this->internalField().objectPath()
244  << exit(FatalError);
245  }
246  }
247 
248  fixedValueFvPatchScalarField::updateCoeffs();
249 }
250 
251 
253 {
255 
256  os.writeEntry("kabs", kabs_);
257  os.writeEntry("max", max_);
258  os.writeEntryIfDifferent<scalar>("kdes", scalar(0), kdes_);
259 
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 namespace Foam
267 {
269  (
271  timeVaryingMassSorptionFvPatchScalarField
272  );
273 }
274 
275 // ************************************************************************* //
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
timeVaryingMassSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: POSIX.C:1063
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
ddtSchemeType
Enumeration defining the available ddt schemes.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
This boundary condition provides a first order fixed-value condition for a given scalar field to mode...
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const dimensionedScalar c
Speed of light in a vacuum.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127