Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
26\*---------------------------------------------------------------------------*/
27
28#include "patchCellsSource.H"
29#include "boundarySourcePatch.H"
30#include "fvMatrices.H"
33// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace fv
38{
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
97 }
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{
168 {
169 return false;
170 }
171
172 return true;
173}
174
175
176// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
if(patchID !=-1)
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....
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
A List with indirect addressing. Like IndirectList but does not store addressing.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
option(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition fvOption.C:51
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition fvOptionIO.C:48
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
Source defined by a boundary condition applied to cells next to patches. This fvOption needs to be us...
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible (momentum, enthalpy, species) equation.
patchCellsSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & psi
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
A special matrix type and solver, designed for finite volume solutions of scalar equations.
auto & name
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299