Loading...
Searching...
No Matches
sensitivityMultiple.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-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "adjointSensitivity.H"
31#include "sensitivityMultiple.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40
43(
44 adjointSensitivity,
45 sensitivityMultiple,
47);
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
51sensitivityMultiple::sensitivityMultiple
52(
53 const fvMesh& mesh,
54 const dictionary& dict,
56)
57:
58 adjointSensitivity(mesh, dict, adjointSolver),
59 sensTypes_(this->dict().get<wordList>("sensitivityTypes")),
60 sens_(sensTypes_.size())
61{
63 {
64 sens_.set
65 (
66 sI,
68 (
69 mesh,
70 this->dict().subDict(sensTypes_[sI]),
72 )
73 );
74 sens_[sI].setSuffix(sensTypes_[sI]);
75 }
76}
77
78
79// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80
82{
84 {
85 forAll(sens_, sI)
86 {
87 sens_[sI].readDict(dict.subDict(sensTypes_[sI]));
88 }
89
90 return true;
91 }
92
93 return false;
94}
95
96
98{
100 {
101 sens_[sI].accumulateIntegrand(dt);
102 }
103}
104
105
107(
108 autoPtr<designVariables>& designVars
109)
110{
112 {
113 sens_[sI].assembleSensitivities(designVars);
114 }
115}
116
117
119(
120 autoPtr<designVariables>& designVars
121)
122{
123 forAll(sens_, sI)
124 {
125 Info<< "Computing sensitivities " << sensTypes_[sI] << endl;
126 sens_[sI].calculateSensitivities(designVars);
128 write(type());
129
130 return derivatives_;
131}
132
133
135{
137 {
138 sens_[sI].clearSensitivities();
139 }
140}
141
142
143void sensitivityMultiple::write(const word& baseName)
144{
145 forAll(sens_, sI)
146 {
147 sens_[sI].write(sensTypes_[sI]);
148 }
149}
150
151
152// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153
154} // End namespace Foam
155
156// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Abstract base class for adjoint-based sensitivities.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
scalarField derivatives_
The sensitivity derivative values.
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
Return a reference to the selected turbulence model.
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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
Calculation of adjoint based sensitivities of multiple types.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
const scalarField & calculateSensitivities(autoPtr< designVariables > &designVars)
Calculates sensitivities at wall surface points.
PtrList< adjointSensitivity > sens_
virtual void write(const word &baseName=word::null)
Write sensitivities to file.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
const fvMesh & mesh() const
Return reference to mesh.
const dictionary & dict() const
Return the construction dictionary.
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
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
runTime write()
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299