speciesSorptionFvPatchScalarField.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) 2022 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 Class
27  Foam::speciesSorptionFvPatchScalarField
28 
29 Group
30  grpGenericBoundaryConditions
31 
32 Description
33  This boundary condition provides a first-order zero-gradient
34  condition for a given scalar field to model time-dependent
35  adsorption-desorption processes.
36 
37  \f[
38  \frac{d c}{d t} = k_{ads} (c_{eq} - c_{abs})
39  \f]
40 
41  where
42  \vartable
43  c_{eq} | Equilibrium concentration
44  c_{abs} | Absorbed at wall
45  k_{ads} | Adsorption rate constant [1/s]
46  \endvartable
47 
48  \f[
49  c_{eq} = c_{max} \frac{k_l \, c_{int}}{1 + k_l \, c_{int}}
50  \f]
51 
52  where
53  \vartable
54  c_{max} | Maximum concentration
55  k_l | Langmuir constant
56  c_{int} | Local cell value concentration
57  \endvartable
58 
59 Usage
60  Example of the boundary condition specification:
61  \verbatim
62  <patchName>
63  {
64  // Mandatory entries
65  type speciesSorption;
66  equilibriumModel <word>;
67  kinematicModel <word>;
68  kabs <scalar>;
69  kl <scalar>;
70  max <scalar>;
71  thickness <PatchFunction1<scalar>>;
72  rhoS <scalar>;
73 
74  // Optional entries
75  dfldp <scalarField>;
76  mass <scalarField>;
77  pName <word>;
78 
79  // Inherited entries
80  ...
81  }
82  \endverbatim
83 
84  where the entries mean:
85  \table
86  Property | Description | Type | Reqd | Deflt
87  type | Type name: speciesSorption | word | yes | -
88  equilibriumModel | Equilibrium model | word | yes | -
89  kinematicModel | Kinematic model | word | yes | -
90  kabs | Adsorption rate constant [1/s] | scalar | yes | -
91  kl | Langmuir constant [1/Pa] | scalar | yes | -
92  max | Maximum concentation at wall [mol/kg] | scalar | yes | -
93  thickness| Solid thickness along the patch <!--
94  --> | PatchFunction1<scalar> | yes | -
95  rhoS | Solid density | scalar | yes | -
96  dfldp | Source on cells next to patch | scalarField | no | Zero
97  mass | Absorbed mass per kg of absorbent [mol/kg] <!--
98  --> | scalarField | no | Zero
99  pName | Name of operand pressure field | word | no | p
100  \endtable
101 
102  Options for the \c equilibriumModel entry:
103  \verbatim
104  Langmuir | Langmuir model
105  \endverbatim
106 
107  Options for the \c kinematicModel entry:
108  \verbatim
109  PseudoFirstOrder | Pseudo first-order model
110  \endverbatim
111 
112  The inherited entries are elaborated in:
113  - \link zeroGradientFvPatchFields.H \endlink
114  - \link PatchFunction1.H \endlink
115 
116 SourceFiles
117  speciesSorptionFvPatchScalarField.C
118 
119 \*---------------------------------------------------------------------------*/
120 
121 #ifndef Foam_speciesSorptionFvPatchScalarField_H
122 #define Foam_speciesSorptionFvPatchScalarField_H
123 
124 #include "boundarySourcePatch.H"
126 #include "PatchFunction1.H"
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
130 namespace Foam
131 {
132 
133 /*---------------------------------------------------------------------------*\
134  Class speciesSorptionFvPatchScalarField Declaration
135 \*---------------------------------------------------------------------------*/
136 
137 class speciesSorptionFvPatchScalarField
138 :
139  public zeroGradientFvPatchField<scalar>,
140  public boundarySourcePatch
141 {
142 public:
143 
144  // Public Enumeration
145 
146  //- Options for the equilibrum model
147  enum equilibriumModelType : char
148  {
149  LANGMUIR = 0
150  };
151 
152  //- Options for the kinematic model
153  enum kineticModelType : char
154  {
155  PseudoFirstOrder = 0
156  };
157 
158  //- Names for equilibriumModelType
159  static const Enum<equilibriumModelType> equilibriumModelTypeNames;
160 
161  //- Names for kineticModelType
162  static const Enum<kineticModelType> kinematicModelTypeNames;
163 
164 
165 private:
166 
167  // Private Data
168 
169  //- Equilibrium model
170  enum equilibriumModelType equilibriumModel_;
171 
172  //- Kinematic model
173  enum kineticModelType kinematicModel_;
174 
175  //- Solid thickness along the patch
176  autoPtr<PatchFunction1<scalar>> thicknessPtr_;
177 
178  //- Adsorption rate constant [1/sec]
179  scalar kabs_;
180 
181  //- Langmuir adsorption constant [1/Pa]
182  scalar kl_;
183 
184  //- Maximum density on patch [mol/kg]
185  scalar max_;
186 
187  //- Solid density
188  scalar rhoS_;
189 
190  //- Name of operand pressure field
191  word pName_;
192 
193  //- Source on cells next to patch [mol/kg/sec]
194  scalarField dfldp_;
195 
196  //- Absorbed mass per kg of absorbent [mol/kg]
197  scalarField mass_;
198 
199 
200  // Private Member Functions
201 
202  //- Calculate the mole fraction fields
203  tmp<scalarField> calcMoleFractions() const;
204 
205  //- Lookup (or create) field for output
206  volScalarField& field(const word&, const dimensionSet&) const;
207 
208 
209 public:
210 
211  //- Runtime type information
212  TypeName("speciesSorption");
213 
214 
215  // Constructors
216 
217  //- Construct from patch and internal field
219  (
220  const fvPatch&,
221  const DimensionedField<scalar, volMesh>&
222  );
223 
224  //- Construct from patch, internal field and dictionary
226  (
227  const fvPatch&,
228  const DimensionedField<scalar, volMesh>&,
229  const dictionary&
230  );
231 
232  //- Construct by mapping given
233  //- speciesSorptionFvPatchScalarField onto a new patch
235  (
237  const fvPatch&,
239  const fvPatchFieldMapper&
240  );
241 
242  //- Construct as copy
244  (
246  );
247 
248  //- Construct as copy setting internal field reference
250  (
253  );
255  //- Return a clone
256  virtual tmp<fvPatchField<scalar>> clone() const
257  {
258  return fvPatchField<scalar>::Clone(*this);
259  }
261  //- Clone with an internal field reference
263  (
265  ) const
266  {
267  return fvPatchField<scalar>::Clone(*this, iF);
268  }
269 
270 
271  // Member Functions
272 
273  // Mapping
274 
275  //- Map (and resize as needed) from self given a mapping object
276  virtual void autoMap
277  (
278  const fvPatchFieldMapper&
279  );
280 
281  //- Reverse map the given fvPatchField onto this fvPatchField
282  virtual void rmap
283  (
284  const fvPatchScalarField&,
285  const labelList&
286  );
287 
288 
289  // Evaluation
290 
291  //- Source of cells next to the patch
292  virtual tmp<scalarField> patchSource() const;
293 
294  //- Access to mass
295  tmp<scalarField> mass() const;
296 
297  //- Update the coefficients associated with the patch field
298  virtual void updateCoeffs();
299 
300 
301  // I-O
302 
303  //- Write
304  virtual void write(Ostream&) const;
305 };
306 
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 } // End namespace Foam
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 #endif
315 
316 // ************************************************************************* //
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
TypeName("speciesSorption")
Runtime type information.
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
speciesSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition provides a first-order zero-gradient condition for a given scalar field to mo...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
Definition: fvPatchField.H:597
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static const Enum< equilibriumModelType > equilibriumModelTypeNames
Names for equilibriumModelType.
tmp< scalarField > mass() const
Access to mass.
static const Enum< kineticModelType > kinematicModelTypeNames
Names for kineticModelType.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual tmp< scalarField > patchSource() const
Source of cells next to the patch.
Namespace for OpenFOAM.