fieldRegularisation.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) 2007-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fieldRegularisation.H"
32 #include "wallFvPatch.H"
33 #include "fvm.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 Foam::fieldRegularisation::fieldRegularisation
45 (
46  fvMesh& mesh,
47  const scalarField& alpha,
48  const topOZones& zones,
49  const dictionary& dict
50 )
51 :
52  mesh_(mesh),
53  dict_(dict),
54  zones_(zones),
55  regularise_(dict.getOrDefault<bool>("regularise", false)),
56  project_(dict.getOrDefault<bool>("project", regularise_)),
57  radius_(regularisationRadius::New(mesh, dict, false)),
58  alpha_(alpha),
59  alphaTilda_
60  (
61  regularise_
62  ? new volScalarField
63  (
64  IOobject
65  (
66  "alphaTilda",
67  mesh_.time().timeName(),
68  mesh_,
69  IOobject::READ_IF_PRESENT,
70  IOobject::AUTO_WRITE
71  ),
72  mesh_,
74  fixedValueFvPatchScalarField::typeName
75  )
76  : nullptr
77  ),
78  sharpenFunction_
79  (
80  project_ ?
82  nullptr
83  ),
84  regularisationPDE_(regularisationPDE::New(mesh, dict, zones)),
85  betaArg_(regularise_ ? alphaTilda_().primitiveField() : alpha_),
86  growFromWalls_(dict.getOrDefault<bool>("growFromWalls", false)),
87  beta_
88  (
89  IOobject
90  (
91  "beta",
92  mesh_.time().timeName(),
93  mesh_,
94  IOobject::READ_IF_PRESENT,
95  IOobject::AUTO_WRITE
96  ),
97  mesh_,
99  zeroGradientFvPatchScalarField::typeName
100  )
101 {
102  DebugInfo
103  << "Regularise " << Switch(regularise_) << nl
104  << "Project " << Switch(project_) << endl;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
112  if (regularise_)
113  {
114  regularise(alpha_, alphaTilda_(), true);
115  }
116 
117  if (project_)
118  {
119  sharpenFunction_->interpolate(betaArg_, beta_.primitiveFieldRef());
120  }
121  else
122  {
123  beta_.primitiveFieldRef() = betaArg_;
124  }
125 
126  beta_.correctBoundaryConditions();
127 }
128 
129 
131 (
132  const scalarField& source,
133  scalarField& result,
134  const bool isTopoField,
135  const regularisationRadius& radius
136 )
137 {
138  regularisationPDE_->
139  regularise(alphaTilda_(), source, result, isTopoField, radius);
140 }
141 
142 
144 (
145  const scalarField& source,
146  scalarField& result,
147  const bool isTopoField
148 )
149 {
150  regularisationPDE_->
151  regularise(alphaTilda_(), source, result, isTopoField, radius_());
152 }
153 
154 
156 {
157  // Add dBeta/dBetaArg
158  if (project_)
159  {
160  sens *= sharpenFunction_->derivative(betaArg_);
161  }
162  // Add part due to regularisation
163  if (regularise_)
164  {
165  // Solve the adjoint to the regularisation equation
166  regularise(sens, sens, false);
167  }
168 
169  // Add volume
170  sens *= mesh_.V();
171 }
172 
173 
174 // ************************************************************************* //
Base class for selecting the regulatisation radius.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool regularise_
Perform regulaisation on alpha before inputing it on beta?
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool project_
Perform the projection (sharpening) step?
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
virtual void updateBeta()
Update the beta field.
Base class for selecting the regulatisation PDE.
const dimensionSet dimless
Dimensionless.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
void regularise(const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius)
Regularise an externally provided radius.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInfo
Report an information message using Foam::Info.
defineTypeNameAndDebug(combustionModel, 0)
void postProcessSens(scalarField &sens)
Update part of fieldRegularisation to the sensitivitiy derivatives.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127