Loading...
Searching...
No Matches
adjointFarFieldPressureFvPatchScalarField.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
33#include "volFields.H"
34#include "surfaceFields.H"
35#include "ATCUaGradU.H"
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
41(
42 const fvPatch& p,
57 const fvPatchFieldMapper& mapper
59:
60 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
62{}
63
64
67(
68 const fvPatch& p,
70 const dictionary& dict
71)
72:
73 fixedValueFvPatchScalarField(p, iF),
74 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
75{
76 this->readValueEntry(dict, IOobjectOption::MUST_READ);
77}
78
79
82(
85)
86:
87 fixedValueFvPatchScalarField(tppsf, iF),
89{}
90
91
92// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93
95{
96 if (updated())
97 {
98 return;
99 }
100
101 // Patch normal and surface
102 const scalarField& magSf = patch().magSf();
103 const vectorField nf(patch().nf());
104
105 // Primal flux
106 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
107
108 // Adjoint flux
109 //const fvsPatchField<scalar>& phiap =
110 // patch().lookupPatchField<surfaceScalarField>("phia");
111
112 // Primal velocity
113 const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
114
115 // Adjoint velocity
116 const fvPatchField<vector>& Uap = boundaryContrPtr_->Uab();
117
118 // Patch-adjacent normal adjoint velocity
119 scalarField Uac_n(Uap.patchInternalField()() & nf);
120
121 // Patch normal adjoint velocity
122 scalarField Uap_n(Uap & nf);
123 //scalarField Uap_n = phiap/magSf;
124
125 // Patch normal velocity Uap_n
126 scalarField phiOverSurf(phip/magSf);
127
128 // Patch deltas
129 const scalarField& delta = patch().deltaCoeffs();
130
131 // snGrad Ua_n
132 scalarField snGradUan(delta*(Uap_n - Uac_n));
133
134 // Momentum diffusion coefficient
135 tmp<scalarField> tmomentumDiffusion =
136 boundaryContrPtr_->momentumDiffusion();
137 scalarField& momentumDiffusion = tmomentumDiffusion.ref();
138
139 // Objective function and other explicit contributions
140 tmp<scalarField> tsource = boundaryContrPtr_->pressureSource();
141 scalarField source = tsource.ref();
142
143 // Contribution from the ATC part (if UaGradU)
144 if (addATCUaGradUTerm())
145 {
146 source += Uap & Up;
147 }
148
149 operator==
150 (
151 // Inlet
152 neg(phip)*(patchInternalField())
153
154 // Outlet
155 + pos(phip)*
156 (
157 Uap_n*phiOverSurf
158 + 2*momentumDiffusion*snGradUan
159 + source
160 )
161 );
162
163 fixedValueFvPatchScalarField::updateCoeffs();
164}
165
166
169{
170 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
171
172 return tmp<Field<scalar>>
173 (
174 new Field<scalar>
175 (
176 pos(phip)*patch().deltaCoeffs()*(*this - patchInternalField())
177 )
178 );
179}
180
181
184(
185 const tmp<scalarField>&
186) const
187{
188 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
189
190 return tmp<Field<scalar>>
191 (
192 new Field<scalar>
193 (
195 )
196 );
197}
198
199
202(
203 const tmp<scalarField>&
204) const
205{
206 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
207
208 return tmp<Field<scalar>>
209 (
210 new Field<scalar>
211 (
212 pos(phip)*(*this)
213 )
214 );
215}
216
217
220{
221 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
222
223 // Act as a zeroGradient pa bc
224 return tmp<Field<scalar>>
225 (
226 new Field<scalar>
227 (
228 -pos(phip)*pTraits<scalar>::one*(this->patch().deltaCoeffs())
229 )
230 );
231}
232
233
236{
237 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
238
239 // Act as a zeroGradient pa bc
240 return tmp<Field<scalar>>
241 (
242 new Field<scalar>
244 pos(phip)*(this->patch().deltaCoeffs()*(*this))
245 )
246 );
247}
248
249
251{
253 os.writeEntry("solverName", adjointSolverName_);
255}
256
257
258// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
259
260void Foam::adjointFarFieldPressureFvPatchScalarField::operator=
261(
262 const UList<scalar>& ul
263)
264{
266 scalarField value(neg(phip)*ul + pos(phip)*(*this));
267
269}
270
271
272void Foam::adjointFarFieldPressureFvPatchScalarField::operator=
273(
274 const fvPatchField<scalar>& ptf
275)
276{
277 check(ptf);
279 scalarField value(neg(phip)*ptf + pos(phip)*(*this));
280
282}
283
284
285void Foam::adjointFarFieldPressureFvPatchScalarField::operator+=
286(
287 const fvPatchField<scalar>& ptf
288)
289{
290 check(ptf);
292 scalarField value(neg(phip)*((*this) + ptf) + pos(phip)*(*this));
293
295}
296
297
298void Foam::adjointFarFieldPressureFvPatchScalarField::operator-=
299(
300 const fvPatchField<scalar>& ptf
301)
302{
303 check(ptf);
305 scalarField value(neg(phip)*((*this) - ptf) + pos(phip)*(*this));
306
308}
309
310
311void Foam::adjointFarFieldPressureFvPatchScalarField::operator*=
312(
313 const fvPatchField<scalar>& ptf
314)
315{
316 if (&patch() != &ptf.patch())
317 {
319 << "Incompatible patches for patch fields"
320 << abort(FatalError);
321 }
322
324 scalarField value(neg(phip)*((*this)*ptf) + pos(phip)*(*this));
325
327}
328
329
330void Foam::adjointFarFieldPressureFvPatchScalarField::operator/=
331(
332 const fvPatchField<scalar>& ptf
333)
334{
335 if (&patch() != &ptf.patch())
336 {
338 << "Incompatible patches for patch fields"
339 << abort(FatalError);
340 }
341
343 scalarField value(neg(phip)*((*this)/ptf) + pos(phip)*(*this));
344
346}
347
348
349void Foam::adjointFarFieldPressureFvPatchScalarField::operator+=
350(
351 const Field<scalar>& tf
352)
353{
355 scalarField value(neg(phip)*((*this) + tf) + pos(phip)*(*this));
356
358}
359
360
361void Foam::adjointFarFieldPressureFvPatchScalarField::operator-=
362(
363 const Field<scalar>& tf
364)
365{
366 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
367 scalarField value(neg(phip)*((*this)-tf) + pos(phip)*(*this));
368
370}
371
372
373void Foam::adjointFarFieldPressureFvPatchScalarField::operator*=
374(
375 const scalarField& tf
376)
377{
378 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
379 scalarField value(neg(phip)*((*this)*tf) + pos(phip)*(*this));
380
382}
383
384
385void Foam::adjointFarFieldPressureFvPatchScalarField::operator/=
386(
387 const scalarField& tf
388)
389{
391 scalarField value(neg(phip)*((*this)/tf) + pos(phip)*(*this));
392
394}
395
396
397void Foam::adjointFarFieldPressureFvPatchScalarField::operator=
398(
399 const scalar t
400)
401{
403 scalarField value(neg(phip)*t + pos(phip)*(*this));
404
406}
407
408
409void Foam::adjointFarFieldPressureFvPatchScalarField::operator+=
410(
411 const scalar t
412)
413{
415 scalarField value(neg(phip)*((*this) + t) + pos(phip)*(*this));
416
418}
419
420
421void Foam::adjointFarFieldPressureFvPatchScalarField::operator-=
422(
423 const scalar t
424)
425{
426 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
427 scalarField value
428 (
429 neg(phip)*((*this)-t)
430 + pos(phip)*(*this)
431 );
432
434}
435
436
437void Foam::adjointFarFieldPressureFvPatchScalarField::operator*=
438(
439 const scalar s
440)
441{
443 scalarField value(neg(phip)*((*this)*s) + pos(phip)*(*this));
444
446}
447
448
449void Foam::adjointFarFieldPressureFvPatchScalarField::operator/=
450(
451 const scalar s
452)
453{
454 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
455 scalarField value(neg(phip)*((*this)/s) + pos(phip)*(*this));
456
459
460
461// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
462
463namespace Foam
464{
466 (
468 adjointFarFieldPressureFvPatchScalarField
469 );
470}
471
472// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void operator=(const Field< Type > &)
Copy assignment.
Definition Field.C:781
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
virtual tmp< Field< scalar > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this patch...
virtual tmp< Field< scalar > > snGrad() const
Return true if this patch field fixes a value.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< Field< scalar > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the evaluation of the gradient of this patchFi...
adjointFarFieldPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
virtual tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
static void check(const int retVal, const char *what)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
fvPatchField< scalar > fvPatchScalarField
dictionary dict
Foam::surfaceFields.