adjointFarFieldTMVar1FvPatchScalarField.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) 2014-2022 PCOpt/NTUA
9  Copyright (C) 2014-2022 FOSS GP
10 -------------------------------------------------------------------------------
11 License
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 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
44 (
45  const fvPatch& p,
47 )
48 :
49  fixedValueFvPatchScalarField(p, iF),
51 {}
52 
53 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
63  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
64  adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
65 {}
66 
67 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  fixedValueFvPatchScalarField(p, iF),
77  adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
78 {
80  (
81  scalarField("value", dict, p.size())
82  );
83 }
84 
85 
88 (
91 )
92 :
93  fixedValueFvPatchScalarField(tppsf, iF),
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (updated())
103  {
104  return;
105  }
106  vectorField nf(patch().nf());
107 
108  tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
109  const scalarField& nuEff = tnuEff();
110  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
111  const scalarField& magSf = patch().magSf();
112 
113  // Patch-adjacent values
114  tmp<scalarField> intf(patchInternalField());
115 
116  // Patch deltas
117  const scalarField& delta = patch().deltaCoeffs();
118 
119  operator==
120  (
121  pos(phip)
122  *(
123  (nuEff*delta*intf)
124  /(phip/magSf + nuEff*delta)
125  )
126  );
128  fixedValueFvPatchScalarField::updateCoeffs();
129 }
130 
131 
134 (
135  const tmp<scalarField>&
136 ) const
137 {
138  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
139 
140  // For fixedValue patches
142 }
143 
144 
147 (
148  const tmp<scalarField>&
149 ) const
150 {
152 
153  // For zeroGradient patches
154  return tmp<Field<scalar>>::New(pos(phip)*(*this));
155 }
156 
157 
159 {
161  writeEntry("value", os);
162  os.writeEntry("solverName", adjointSolverName_);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
169 (
171  adjointFarFieldTMVar1FvPatchScalarField
172 );
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 } // End namespace Foam
177 
178 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
scalar delta
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
A traits class, which is primarily used for primitives.
Definition: pTraits.H:50
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
adjointFarFieldTMVar1FvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensionedScalar pos(const dimensionedScalar &ds)
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Engine to manage contributions of the objective functions to the adjoint boundary conditions...
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
static const word null
An empty word.
Definition: word.H:84
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
OBJstream os(runTime.globalPath()/outputName)
Base class for solution control classes.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
virtual tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
word adjointSolverName_
adjointSolver name corresponding to field
Namespace for OpenFOAM.