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, false),
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 (dict.found("value"))
89  {
90  fvPatchScalarField::operator=
91  (
92  scalarField("value", dict, p.size())
93  );
94  }
95  else
96  {
98  }
99 }
100 
101 
104 (
106  const fvPatch& p,
108  const fvPatchFieldMapper& mapper
109 )
110 :
111  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
112  kabs_(ptf.kabs_),
113  max_(ptf.max_),
114  kdes_(ptf.kdes_)
115 {}
116 
117 
120 (
122 )
123 :
124  fixedValueFvPatchScalarField(ptf),
125  kabs_(ptf.kabs_),
126  max_(ptf.max_),
127  kdes_(ptf.kdes_)
128 {}
129 
130 
133 (
136 )
137 :
138  fixedValueFvPatchScalarField(ptf, iF),
139  kabs_(ptf.kabs_),
140  max_(ptf.max_),
141  kdes_(ptf.kdes_)
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
148 (
149  const fvPatchFieldMapper& m
150 )
151 {
152  fixedValueFvPatchScalarField::autoMap(m);
153 }
154 
155 
157 (
158  const fvPatchScalarField& ptf,
159  const labelList& addr
160 )
161 {
162  fixedValueFvPatchScalarField::rmap(ptf, addr);
163 }
164 
165 
168 {
169  auto tsource = tmp<scalarField>::New(patch().size(), Zero);
170  auto& source = tsource.ref();
171 
172  const scalarField cp(*this);
173  const scalarField w(max(1 - cp/max_, scalar(0)));
174 
175  source = -kabs_*w*max(patchInternalField() - cp, scalar(0));
176 
177  source += kdes_*max(cp - patchInternalField(), scalar(0));
178 
179  return tsource;
180 }
181 
182 
184 {
185  if (updated())
186  {
187  return;
188  }
189 
190  const label patchi = patch().index();
191 
192  const scalar dt = db().time().deltaTValue();
193 
194  const auto& fld =
195  db().lookupObject<volScalarField>(this->internalField().name());
196  const volScalarField& fld0 = fld.oldTime();
197 
198  // Lookup d/dt scheme from database
199  const word ddtSchemeName(fld.mesh().ddtScheme(fld.name()));
200  const ddtSchemeType& ddtScheme = ddtSchemeTypeNames_[ddtSchemeName];
201 
202  const scalarField cp(*this);
203  const scalarField w(max(1 - cp/max_, scalar(0)));
204 
205  tmp<scalarField> dfldp =
206  kabs_
207  *w
208  *max(patchInternalField() - cp, scalar(0))
209  *dt;
210 
211  dfldp.ref() -=
212  kdes_
213  *max(cp - patchInternalField(), scalar(0))
214  *dt;
215 
216  switch (ddtScheme)
217  {
218  case tsEuler:
219  case tsCrankNicolson:
220  {
221  operator==(fld0.boundaryField()[patchi] + dfldp);
222 
223  break;
224  }
225  case tsBackward:
226  {
227  const scalar dt0 = db().time().deltaT0Value();
228 
229  const scalar c = scalar(1) + dt/(dt + dt0);
230  const scalar c00 = dt*dt/(dt0*(dt + dt0));
231  const scalar c0 = c + c00;
232 
233  operator==
234  (
235  (
236  c0*fld0.boundaryField()[patchi]
237  - c00*fld0.oldTime().boundaryField()[patchi]
238  + dfldp
239  )/c
240  );
241 
242  break;
243  }
244  default:
245  {
247  << ddtSchemeName << nl
248  << " on patch " << this->patch().name()
249  << " of field " << this->internalField().name()
250  << " in file " << this->internalField().objectPath()
251  << exit(FatalError);
252  }
253  }
254 
255  fixedValueFvPatchScalarField::updateCoeffs();
256 }
257 
258 
260 {
262 
263  os.writeEntry("kabs", kabs_);
264  os.writeEntry("max", max_);
265  os.writeEntryIfDifferent<scalar>("kdes", scalar(0), kdes_);
266 
267  writeEntry("value", os);
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 namespace Foam
274 {
276  (
279  );
280 }
281 
282 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
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:248
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
backwardDdtScheme.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
fvPatchFieldMapper.H
Foam::timeVaryingMassSorptionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: timeVaryingMassSorptionFvPatchScalarField.C:259
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::timeVaryingMassSorptionFvPatchScalarField::ddtSchemeType
ddtSchemeType
Enumeration defining the available ddt schemes.
Definition: timeVaryingMassSorptionFvPatchScalarField.H:167
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
EulerDdtScheme.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::timeVaryingMassSorptionFvPatchScalarField::timeVaryingMassSorptionFvPatchScalarField
timeVaryingMassSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: timeVaryingMassSorptionFvPatchScalarField.C:63
Foam::timeVaryingMassSorptionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: timeVaryingMassSorptionFvPatchScalarField.C:148
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
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))
Definition: gmvOutputLagrangian.H:23
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
CrankNicolsonDdtScheme.H
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::timeVaryingMassSorptionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: timeVaryingMassSorptionFvPatchScalarField.C:183
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::cp
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: MSwindows.C:802
Foam::nl
constexpr char nl
Definition: Ostream.H:404
timeVaryingMassSorptionFvPatchScalarField.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::timeVaryingMassSorptionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: timeVaryingMassSorptionFvPatchScalarField.C:157
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::timeVaryingMassSorptionFvPatchScalarField
This boundary condition provides a first order fixed-value condition for a given scalar field to mode...
Definition: timeVaryingMassSorptionFvPatchScalarField.H:158
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::timeVaryingMassSorptionFvPatchScalarField::source
tmp< scalarField > source() const
Return source rate.
Definition: timeVaryingMassSorptionFvPatchScalarField.C:167
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54