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 -------------------------------------------------------------------------------
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 
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 #include "ATCUaGradU.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  fixedValueFvPatchScalarField(p, iF),
47  adjointScalarBoundaryCondition(p, iF, word::null)
48 {}
49 
50 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
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 {
77  (
78  scalarField("value", dict, p.size())
79  );
80 }
81 
82 
85 (
88 )
89 :
90  fixedValueFvPatchScalarField(tppsf, iF),
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  if (updated())
100  {
101  return;
102  }
103 
104  // Patch normal and surface
105  const scalarField& magSf = patch().magSf();
106  const vectorField nf(patch().nf());
107 
108  // Primal flux
109  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
110 
111  // Adjoint flux
112  //const fvsPatchField<scalar>& phiap =
113  // patch().lookupPatchField<surfaceScalarField, scalar>("phia");
114 
115  // Primal velocity
116  const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
117 
118  // Adjoint velocity
119  const fvPatchField<vector>& Uap = boundaryContrPtr_->Uab();
120 
121  // Patch-adjacent normal adjoint velocity
122  scalarField Uac_n(Uap.patchInternalField()() & nf);
123 
124  // Patch normal adjoint velocity
125  scalarField Uap_n(Uap & nf);
126  //scalarField Uap_n = phiap/magSf;
127 
128  // Patch normal velocity Uap_n
129  scalarField phiOverSurf(phip/magSf);
130 
131  // Patch deltas
132  const scalarField& delta = patch().deltaCoeffs();
133 
134  // snGrad Ua_n
135  scalarField snGradUan(delta*(Uap_n - Uac_n));
136 
137  // Momentum diffusion coefficient
138  tmp<scalarField> tmomentumDiffusion =
139  boundaryContrPtr_->momentumDiffusion();
140  scalarField& momentumDiffusion = tmomentumDiffusion.ref();
141 
142  // Objective function and other explicit contributions
143  tmp<scalarField> tsource = boundaryContrPtr_->pressureSource();
144  scalarField source = tsource.ref();
145 
146  // Contribution from the ATC part (if UaGradU)
147  if (addATCUaGradUTerm())
148  {
149  source += Uap & Up;
150  }
151 
152  operator==
153  (
154  // Inlet
155  neg(phip)*(patchInternalField())
156 
157  // Outlet
158  + pos(phip)*
159  (
160  Uap_n*phiOverSurf
161  + 2*momentumDiffusion*snGradUan
162  + source
163  )
164  );
165 
166  fixedValueFvPatchScalarField::updateCoeffs();
167 }
168 
169 
172 {
173  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
174 
175  return tmp<Field<scalar>>
176  (
177  new Field<scalar>
178  (
179  pos(phip)*patch().deltaCoeffs()*(*this - patchInternalField())
180  )
181  );
182 }
183 
184 
187 (
188  const tmp<scalarField>&
189 ) const
190 {
191  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
192 
193  return tmp<Field<scalar>>
194  (
195  new Field<scalar>
196  (
198  )
199  );
200 }
201 
202 
205 (
206  const tmp<scalarField>&
207 ) const
208 {
209  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
210 
211  return tmp<Field<scalar>>
212  (
213  new Field<scalar>
214  (
215  pos(phip)*(*this)
216  )
217  );
218 }
219 
220 
223 {
224  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
225 
226  // Act as a zeroGradient pa bc
227  return tmp<Field<scalar>>
228  (
229  new Field<scalar>
230  (
231  -pos(phip)*pTraits<scalar>::one*(this->patch().deltaCoeffs())
232  )
233  );
234 }
235 
236 
239 {
240  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
241 
242  // Act as a zeroGradient pa bc
243  return tmp<Field<scalar>>
244  (
245  new Field<scalar>
246  (
247  pos(phip)*(this->patch().deltaCoeffs()*(*this))
248  )
249  );
250 }
251 
252 
254 {
256  writeEntry("value", os);
257  os.writeEntry("solverName", adjointSolverName_);
258 }
259 
260 
261 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
262 
263 void Foam::adjointFarFieldPressureFvPatchScalarField::operator=
264 (
265  const UList<scalar>& ul
266 )
267 {
268  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
269  scalarField value(neg(phip)*ul + pos(phip)*(*this));
270 
272 }
273 
274 
275 void Foam::adjointFarFieldPressureFvPatchScalarField::operator=
276 (
277  const fvPatchField<scalar>& ptf
278 )
279 {
280  check(ptf);
281  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
282  scalarField value(neg(phip)*ptf + pos(phip)*(*this));
283 
285 }
286 
287 
288 void Foam::adjointFarFieldPressureFvPatchScalarField::operator+=
289 (
290  const fvPatchField<scalar>& ptf
291 )
292 {
293  check(ptf);
294  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
295  scalarField value(neg(phip)*((*this) + ptf) + pos(phip)*(*this));
296 
298 }
299 
300 
301 void Foam::adjointFarFieldPressureFvPatchScalarField::operator-=
302 (
303  const fvPatchField<scalar>& ptf
304 )
305 {
306  check(ptf);
307  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
308  scalarField value(neg(phip)*((*this) - ptf) + pos(phip)*(*this));
309 
311 }
312 
313 
314 void Foam::adjointFarFieldPressureFvPatchScalarField::operator*=
315 (
316  const fvPatchField<scalar>& ptf
317 )
318 {
319  if (&patch() != &ptf.patch())
320  {
322  << "Incompatible patches for patch fields"
323  << abort(FatalError);
324  }
325 
326  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
327  scalarField value(neg(phip)*((*this)*ptf) + pos(phip)*(*this));
328 
330 }
331 
332 
333 void Foam::adjointFarFieldPressureFvPatchScalarField::operator/=
334 (
335  const fvPatchField<scalar>& ptf
336 )
337 {
338  if (&patch() != &ptf.patch())
339  {
341  << "Incompatible patches for patch fields"
342  << abort(FatalError);
343  }
344 
345  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
346  scalarField value(neg(phip)*((*this)/ptf) + pos(phip)*(*this));
347 
349 }
350 
351 
352 void Foam::adjointFarFieldPressureFvPatchScalarField::operator+=
353 (
354  const Field<scalar>& tf
355 )
356 {
357  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
358  scalarField value(neg(phip)*((*this) + tf) + pos(phip)*(*this));
359 
361 }
362 
363 
364 void Foam::adjointFarFieldPressureFvPatchScalarField::operator-=
365 (
366  const Field<scalar>& tf
367 )
368 {
369  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
370  scalarField value(neg(phip)*((*this)-tf) + pos(phip)*(*this));
371 
373 }
374 
375 
376 void Foam::adjointFarFieldPressureFvPatchScalarField::operator*=
377 (
378  const scalarField& tf
379 )
380 {
381  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
382  scalarField value(neg(phip)*((*this)*tf) + pos(phip)*(*this));
383 
385 }
386 
387 
388 void Foam::adjointFarFieldPressureFvPatchScalarField::operator/=
389 (
390  const scalarField& tf
391 )
392 {
393  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
394  scalarField value(neg(phip)*((*this)/tf) + pos(phip)*(*this));
395 
397 }
398 
399 
400 void Foam::adjointFarFieldPressureFvPatchScalarField::operator=
401 (
402  const scalar t
403 )
404 {
405  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
406  scalarField value(neg(phip)*t + pos(phip)*(*this));
407 
409 }
410 
411 
412 void Foam::adjointFarFieldPressureFvPatchScalarField::operator+=
413 (
414  const scalar t
415 )
416 {
417  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
418  scalarField value(neg(phip)*((*this) + t) + pos(phip)*(*this));
419 
421 }
422 
423 
424 void Foam::adjointFarFieldPressureFvPatchScalarField::operator-=
425 (
426  const scalar t
427 )
428 {
429  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
430  scalarField value
431  (
432  neg(phip)*((*this)-t)
433  + pos(phip)*(*this)
434  );
435 
437 }
438 
439 
440 void Foam::adjointFarFieldPressureFvPatchScalarField::operator*=
441 (
442  const scalar s
443 )
444 {
445  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
446  scalarField value(neg(phip)*((*this)*s) + pos(phip)*(*this));
447 
449 }
450 
451 
452 void Foam::adjointFarFieldPressureFvPatchScalarField::operator/=
453 (
454  const scalar s
455 )
456 {
457  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
458  scalarField value(neg(phip)*((*this)/s) + pos(phip)*(*this));
459 
461 }
462 
463 
464 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
465 
466 namespace Foam
467 {
469  (
472  );
473 }
474 
475 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
s
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))
Definition: gmvOutputSpray.H:25
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
surfaceFields.H
Foam::surfaceFields.
Foam::adjointBoundaryCondition::boundaryContrPtr_
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Definition: adjointBoundaryCondition.H:73
fvPatchFieldMapper.H
Foam::adjointFarFieldPressureFvPatchScalarField::gradientInternalCoeffs
virtual tmp< Field< scalar > > gradientInternalCoeffs() const
Definition: adjointFarFieldPressureFvPatchScalarField.C:222
Foam::adjointScalarBoundaryCondition
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
Definition: adjointBoundaryConditionsFwd.H:41
Foam::adjointFarFieldPressureFvPatchScalarField::valueBoundaryCoeffs
virtual tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: adjointFarFieldPressureFvPatchScalarField.C:205
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
ATCUaGradU.H
Foam::Field< scalar >
Foam::adjointFarFieldPressureFvPatchScalarField::gradientBoundaryCoeffs
virtual tmp< Field< scalar > > gradientBoundaryCoeffs() const
Definition: adjointFarFieldPressureFvPatchScalarField.C:238
Foam::adjointBoundaryCondition::adjointSolverName_
word adjointSolverName_
adjointSolver name corresponding to field
Definition: adjointBoundaryCondition.H:65
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
Foam::adjointFarFieldPressureFvPatchScalarField::adjointFarFieldPressureFvPatchScalarField
adjointFarFieldPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: adjointFarFieldPressureFvPatchScalarField.C:41
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::adjointBoundaryCondition::addATCUaGradUTerm
bool addATCUaGradUTerm()
Whether to add the extra term from the UaGradU formulation.
Definition: adjointBoundaryCondition.C:144
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::Field::operator=
void operator=(const Field< Type > &)
Copy assignment.
Definition: Field.C:635
Foam::adjointFarFieldPressureFvPatchScalarField::valueInternalCoeffs
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
Definition: adjointFarFieldPressureFvPatchScalarField.C:187
Foam::adjointFarFieldPressureFvPatchScalarField::snGrad
virtual tmp< Field< scalar > > snGrad() const
Return true if this patch field fixes a value.
Definition: adjointFarFieldPressureFvPatchScalarField.C:171
Foam::adjointFarFieldPressureFvPatchScalarField
Definition: adjointFarFieldPressureFvPatchScalarField.H:54
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::adjointFarFieldPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: adjointFarFieldPressureFvPatchScalarField.C:253
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::adjointFarFieldPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: adjointFarFieldPressureFvPatchScalarField.C:97
adjointFarFieldPressureFvPatchScalarField.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177