timeVaryingMassSorptionFvPatchScalarField.H
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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::timeVaryingMassSorptionFvPatchScalarField
28
29Group
30 grpGenericBoundaryConditions
31
32Description
33 This boundary condition provides a first order fixed-value condition
34 for a given scalar field to model time-dependent adsorption-desoprtion
35 processes to be used with the \c interfaceOxideRate mass model
36
37 \f[
38 \frac{d c}{d t} =
39 k_{abs} w (c_{int} - c_{p_{w}}) + k_{des} (c_{p_{w}} - c_{int})
40 \f]
41
42 \f[
43 w = \max(1 - c_{p_{w}}/max, 0)
44 \f]
45
46 where
47 \vartable
48 c_{int} | Concentration at cell
49 c_{p_{w}} | Concentration at wall
50 k_{abs} | Adsorption rate constant [1/s]
51 k_{des} | Desorption rate constant [1/s]
52 w | Weight function
53 max | Max concentration at wall
54 \endvartable
55
56Usage
57 Example of the boundary condition specification:
58 \verbatim
59 <patchName>
60 {
61 // Mandatory entries
62 type timeVaryingMassSorption;
63 kbas <scalar>;
64 max <scalar>;
65
66 // Optional entries
67 kdes <scalar>;
68
69 // Inherited entries
70 ...
71 }
72 \endverbatim
73
74 where the entries mean:
75 \table
76 Property | Description | Type | Reqd | Deflt
77 type | Type name: timeVaryingAdsorption | word | yes | -
78 kbas | Adsorption rate constant | scalar | yes | -
79 max | Maximum concentation at wall | scalar | yes | -
80 kdes | Desorption rate constant | scalar | no | 0
81 \endtable
82
83 The inherited entries are elaborated in:
84 - \link fixedValueFvPatchFields.H \endlink
85
86SourceFiles
87 timeVaryingMassSorptionFvPatchScalarField.C
88
89\*---------------------------------------------------------------------------*/
90
91#ifndef timeVaryingMassSorptionFvPatchScalarField_H
92#define timeVaryingMassSorptionFvPatchScalarField_H
93
95
96// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98namespace Foam
99{
100
101/*---------------------------------------------------------------------------*\
102 Class timeVaryingMassSorptionFvPatchScalarField Declaration
103\*---------------------------------------------------------------------------*/
104
105class timeVaryingMassSorptionFvPatchScalarField
106:
107 public fixedValueFvPatchScalarField
108{
109public:
110
111 // Public Enumeration
112
113 //- Enumeration defining the available ddt schemes
114 enum ddtSchemeType
115 {
116 tsEuler,
119 };
120
121
122private:
123
124 // Private Data
125
126 //- Time scheme type names
127 static const Enum<ddtSchemeType> ddtSchemeTypeNames_;
128
129 //- Adsorption rate constant
130 scalar kabs_;
131
132 //- Maximum level of adsorption of a given substance on patch
133 scalar max_;
134
135 //- Desorption rate constant
136 scalar kdes_;
137
138
139public:
140
141 //- Runtime type information
142 TypeName("timeVaryingMassSorption");
143
144
145 // Constructors
146
147 //- Construct from patch and internal field
149 (
150 const fvPatch&,
151 const DimensionedField<scalar, volMesh>&
152 );
153
154 //- Construct from patch, internal field and dictionary
156 (
157 const fvPatch&,
159 const dictionary&
160 );
161
162 //- Construct by mapping given
163 //- timeVaryingMassSorptionFvPatchScalarField onto a new patch
165 (
167 const fvPatch&,
170 );
171
172 //- Construct as copy
174 (
176 );
177
178 //- Construct and return a clone
179 virtual tmp<fvPatchScalarField> clone() const
180 {
182 (
184 );
185 }
186
187 //- Construct as copy setting internal field reference
189 (
192 );
193
194 //- Construct and return a clone setting internal field reference
196 (
198 ) const
199 {
201 (
203 );
204 }
205
206
207 // Member Functions
208
209 // Mapping
210
211 //- Map (and resize as needed) from self given a mapping object
212 virtual void autoMap
213 (
214 const fvPatchFieldMapper&
215 );
216
217 //- Reverse map the given fvPatchField onto this fvPatchField
218 virtual void rmap
219 (
220 const fvPatchScalarField&,
221 const labelList&
222 );
223
224 // Help
225
226 //- Return source rate
227 tmp<scalarField> source() const;
228
229
230 // Evaluation
231
232 //- Update the coefficients associated with the patch field
233 virtual void updateCoeffs();
234
235
236 //- Write
237 virtual void write(Ostream&) const;
238};
239
240
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243} // End namespace Foam
244
245// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246
247#endif
249// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a first order fixed-value condition for a given scalar field to mode...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
Construct and return a clone setting internal field reference.
TypeName("timeVaryingMassSorption")
Runtime type information.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
timeVaryingMassSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: tmp.H:65
Namespace for OpenFOAM.
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73