Loading...
Searching...
No Matches
Helmholtz.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-2023 PCOpt/NTUA
9 Copyright (C) 2021-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
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
28Class
29 Foam::Helmholtz
30
31Description
32 A regulatisation PDE based on a Helmholtz filter
33
34
35SourceFiles
36 Helmholtz.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef Helmholtz_H
41#define Helmholtz_H
42
43#include "regularisationPDE.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
50/*---------------------------------------------------------------------------*\
51 Class Helmholtz Declaration
52\*---------------------------------------------------------------------------*/
53
54class Helmholtz
55:
56 public regularisationPDE
57{
58
59private:
60
61 // Private Member Functions
62
63 //- No copy construct
64 Helmholtz(const Helmholtz&) = delete;
65
66 //- No copy assignment
67 void operator=(const Helmholtz&) = delete;
68
69
70protected:
71
72 // Protected Data
73
74 //- Solve the regularisationPDE only on a subset mesh made of the
75 //- active cell zones
78 //- Fixed value at wall boundaries. Defaults to 1
79 scalar wallValue_;
80
81
82 // Protected Member Functions
83
84 //- Solve the regulatisation equation
85 // The mesh used is the one of aTilda
86 void solveEqn
87 (
88 const volScalarField& aTilda,
89 const scalarField& source,
90 scalarField& result,
91 const bool isTopoField,
92 const regularisationRadius& radius,
93 const scalar minSetValue = Zero,
94 const bool fixATildaValues = true
95 );
96
97
98public:
99
100 //- Runtime type information
101 TypeName("Helmholtz");
102
103
104 // Constructors
105
106 //- Construct from components
107 Helmholtz
108 (
109 const fvMesh& mesh,
110 const dictionary& dict,
111 const topOZones& zones
112 );
113
114
115 //- Destructor
116 virtual ~Helmholtz() = default;
117
118
119 // Member Functions
120
121 //- Regularise field (or sensitivities) using a Laplace PDE
122 virtual void regularise
123 (
124 const volScalarField& aTilda,
125 const scalarField& source,
126 scalarField& result,
127 const bool isTopoField,
128 const regularisationRadius& radius,
129 const scalar minSetValue = Zero,
130 const bool fixATildaValues = true
131 );
132};
133
134
135// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136
137} // End namespace Foam
138
139// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140
141#endif
142
143// ************************************************************************* //
scalar wallValue_
Fixed value at wall boundaries. Defaults to 1.
Definition Helmholtz.H:82
virtual ~Helmholtz()=default
Destructor.
TypeName("Helmholtz")
Runtime type information.
void solveEqn(const volScalarField &aTilda, const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius, const scalar minSetValue=Zero, const bool fixATildaValues=true)
Solve the regulatisation equation.
Definition Helmholtz.C:44
virtual void regularise(const volScalarField &aTilda, const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius, const scalar minSetValue=Zero, const bool fixATildaValues=true)
Regularise field (or sensitivities) using a Laplace PDE.
Definition Helmholtz.C:166
bool solveOnActiveCells_
Solve the regularisationPDE only on a subset mesh made of the active cell zones.
Definition Helmholtz.H:77
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for selecting the regulatisation radius.
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68