fvOptionAdjointListTemplates.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class Type>
34 (
36 )
37 {
38  const word& fieldName = fld.name();
39 
40  forAll(*this, i)
41  {
42  optionAdjoint& source = this->operator[](i);
43 
44  label fieldI = source.applyToField(fieldName);
45 
46  if (fieldI != -1)
47  {
48  source.setApplied(fieldI);
49 
50  if (source.isActive())
51  {
52  if (debug)
53  {
54  Info<< "Correcting source " << source.name()
55  << " for field " << fieldName << endl;
56  }
57 
58  source.correct(fld);
59  }
60  }
61  }
62 }
63 
64 
65 template<class Type>
66 Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionAdjointList::operator()
67 (
69 )
70 {
71  return this->operator()(fld, fld.name());
72 }
73 
74 
75 template<class Type>
76 Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionAdjointList::operator()
77 (
79  const word& fieldName
80 )
81 {
82  checkApplied();
83 
84  const dimensionSet ds = fld.dimensions()/dimTime*dimVolume;
85 
86  tmp<fvMatrix<Type>> tmtx(new fvMatrix<Type>(fld, ds));
87  fvMatrix<Type>& mtx = tmtx.ref();
88 
89  forAll(*this, i)
90  {
91  optionAdjoint& source = this->operator[](i);
92 
93  label fieldI = source.applyToField(fieldName);
94 
95  if (fieldI != -1)
96  {
97  source.setApplied(fieldI);
98 
99  if (source.isActive())
100  {
101  if (debug)
102  {
103  Info<< "Applying source " << source.name() << " to field "
104  << fieldName << endl;
105  }
106 
107  source.addSup(mtx, fieldI);
108  }
109  }
110  }
111 
112  return tmtx;
113 }
114 
115 
116 template<class Type>
117 Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionAdjointList::operator()
118 (
119  const volScalarField& rho,
121 )
122 {
123  return this->operator()(rho, fld, fld.name());
124 }
125 
126 
127 template<class Type>
128 Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionAdjointList::operator()
129 (
130  const volScalarField& rho,
132  const word& fieldName
133 )
134 {
135  checkApplied();
136 
137  const dimensionSet ds = rho.dimensions()*fld.dimensions()/dimTime*dimVolume;
138 
139  tmp<fvMatrix<Type>> tmtx(new fvMatrix<Type>(fld, ds));
140  fvMatrix<Type>& mtx = tmtx.ref();
141 
142  forAll(*this, i)
143  {
144  optionAdjoint& source = this->operator[](i);
145 
146  label fieldI = source.applyToField(fieldName);
147 
148  if (fieldI != -1)
149  {
150  source.setApplied(fieldI);
151 
152  if (source.isActive())
153  {
154  if (debug)
155  {
156  Info<< "Applying source " << source.name() << " to field "
157  << fieldName << endl;
158  }
159 
160  source.addSup(rho, mtx, fieldI);
161  }
162  }
163  }
164 
165  return tmtx;
166 }
167 
168 
169 template<class Type>
170 Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionAdjointList::operator()
171 (
172  const volScalarField& alpha,
173  const volScalarField& rho,
175 )
176 {
177  return this->operator()(alpha, rho, fld, fld.name());
178 }
179 
180 
181 template<class Type>
182 Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionAdjointList::operator()
183 (
184  const volScalarField& alpha,
185  const volScalarField& rho,
187  const word& fieldName
188 )
189 {
190  checkApplied();
191 
192  const dimensionSet ds =
193  alpha.dimensions()*rho.dimensions()*fld.dimensions()/dimTime*dimVolume;
194 
195  tmp<fvMatrix<Type>> tmtx(new fvMatrix<Type>(fld, ds));
196  fvMatrix<Type>& mtx = tmtx.ref();
197 
198  forAll(*this, i)
199  {
200  optionAdjoint& source = this->operator[](i);
201 
202  label fieldI = source.applyToField(fieldName);
203 
204  if (fieldI != -1)
205  {
206  source.setApplied(fieldI);
207 
208  if (source.isActive())
209  {
210  if (debug)
211  {
212  Info<< "Applying source " << source.name() << " to field "
213  << fieldName << endl;
214  }
215 
216  source.addSup(alpha, rho, mtx, fieldI);
217  }
218  }
219  }
220 
221  return tmtx;
222 }
223 
224 
225 template<class Type>
227 {
228  checkApplied();
229 
230  forAll(*this, i)
231  {
232  optionAdjoint& source = this->operator[](i);
233 
234  label fieldI = source.applyToField(eqn.psi().name());
235 
236  if (fieldI != -1)
237  {
238  source.setApplied(fieldI);
239 
240  if (source.isActive())
241  {
242  if (debug)
243  {
244  Info<< "Applying constraint " << source.name()
245  << " to field " << eqn.psi().name() << endl;
246  }
247 
248  source.constrain(eqn, fieldI);
249  }
250  }
251  }
252 }
253 
254 
255 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fv::option::correct
virtual void correct(volScalarField &field)
Definition: fvOption.C:307
Foam::fv::option::setApplied
void setApplied(const label fieldi)
Set the applied flag to true for field index fieldi.
Definition: fvOptionI.H:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::fv::option::name
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:30
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::fv::optionAdjointList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &fld)
Correct.
Definition: fvOptionAdjointListTemplates.C:34
rho
rho
Definition: readInitialConditions.H:96
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:285
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::fv::option::applyToField
virtual label applyToField(const word &fieldName) const
Return index of field name if found in fieldNames list.
Definition: fvOption.C:114
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::fv::option::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Definition: fvOption.C:135
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fv::optionAdjointList::checkApplied
void checkApplied() const
Check that all sources have been applied.
Definition: fvOptionAdjointList.C:77
Foam::fv::optionAdjoint
Similar to fv::option but with additional functionality to contribute to the sensitivity deriavtives.
Definition: fvOptionAdjoint.H:58
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::fv::option::constrain
virtual void constrain(fvMatrix< scalar > &eqn, const label fieldi)
Definition: fvOption.C:279
Foam::fv::option::isActive
virtual bool isActive()
Is the source active?
Definition: fvOption.C:108
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fv::optionAdjointList::constrain
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Definition: fvOptionAdjointListTemplates.C:226
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:390
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< Type, fvPatchField, volMesh >