patchCellsSource.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) 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 \*---------------------------------------------------------------------------*/
27 
28 #include "patchCellsSource.H"
29 #include "boundarySourcePatch.H"
30 #include "fvMatrices.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39  defineTypeNameAndDebug(patchCellsSource, 0);
40  addToRunTimeSelectionTable(option, patchCellsSource, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const word& sourceName,
50  const word& modelType,
51  const dictionary& dict,
52  const fvMesh& mesh
53 )
54 :
55  fv::option(sourceName, modelType, dict, mesh),
56  curTimeIndex_(-1),
57  isEnergySource_(false)
58 {
60 
61  label nFields = 0;
62 
63  if
64  (
66  && fieldNames_[0] != "none"
67  )
68  {
69  ++nFields;
70  }
71 
72  if
73  (
75  && fieldNames_[0] != "none"
76  )
77  {
78  isEnergySource_ = true;
79  ++nFields;
80  }
81 
82  if
83  (
84  coeffs_.readIfPresent("species", fieldNames_[0])
85  && fieldNames_[0] != "none"
86  )
87  {
88  ++nFields;
89  }
90 
91  if (nFields != 1)
92  {
94  << "Must be specified for one field (U | he | species), but "
95  << nFields << " fields were specified!" << endl
96  << exit(FatalIOError);
97  }
98 
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 (
107  const volScalarField& rho,
108  fvMatrix<scalar>& eqn,
109  const label fieldi
110 )
111 {
112  if (curTimeIndex_ == mesh_.time().timeIndex())
113  {
114  return;
115  }
116  curTimeIndex_ = mesh_.time().timeIndex();
117 
118  if (debug)
119  {
120  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
121  }
122 
124  (
125  name_ + eqn.psi().name() + "Su",
126  mesh_,
128  );
129 
130  // If source applied to he, we need to loop over T for BC's
131  const volScalarField& psi =
132  (
133  isEnergySource_
134  ? mesh_.lookupObject<volScalarField>("T")
135  : mesh_.lookupObject<volScalarField>(eqn.psi().name())
136  );
137  const auto& psib = psi.boundaryField();
138 
139  forAll(psib, patchi)
140  {
141  const auto* bpatchPtr = isA<boundarySourcePatch>(psib[patchi]);
142 
143  if (bpatchPtr)
144  {
145  tmp<scalarField> tsbnd = bpatchPtr->patchSource();
146  const auto& sbnd = tsbnd.cref();
147 
148  UIndirectList<scalar>
149  (
150  tsu.ref().field(),
151  mesh_.boundary()[patchi].faceCells()
152  ) = sbnd;
153  }
154  }
155 
156  if (debug)
157  {
158  Info<< "Field source rate min/max : " << gMinMax(tsu()) << endl;
159  }
160 
161  eqn += tsu;
162 }
163 
164 
166 {
167  if (!fv::option::read(dict))
168  {
169  return false;
170  }
171 
172  return true;
173 }
174 
175 
176 // ************************************************************************* //
patchCellsSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
dictionary dict
virtual bool read(const dictionary &dict)
Read source dictionary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:48
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
const dimensionSet & dimensions() const noexcept
Definition: fvMatrix.H:528
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field...
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const volScalarField & psi
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible (momentum, enthalpy, species) equation. ...
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127