adjointFarFieldNuaTildaFvPatchScalarField.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 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
45 (
46  const fvPatch& p,
48 )
49 :
50  fixedValueFvPatchScalarField(p, iF),
52 {}
53 
54 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  fixedValueFvPatchScalarField(p, iF),
78  adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
79 {
81  (
82  scalarField("value", dict, p.size())
83  );
84 }
85 
86 
89 (
92 )
93 :
94  fixedValueFvPatchScalarField(tppsf, iF),
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  if (updated())
104  {
105  return;
106  }
107  vectorField nf(patch().nf());
108 
109  const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
110  tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
111  const scalarField& nuEff = tnuEff();
112  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
113 
114  // Patch-adjacent nuaTilda nuaTildaNei
115  tmp<scalarField> tnuaTildaNei(patchInternalField());
116  const scalarField& nuaTildaNei = tnuaTildaNei();
117 
118  // Patch deltas
119  const scalarField& delta = patch().deltaCoeffs();
120 
121  operator==
122  (
123  pos(phip)
124  *(
125  (nuEff*delta*nuaTildaNei)
126  /((Ub & nf) + nuEff*delta)
127  )
128  );
129 
130  fixedValueFvPatchScalarField::updateCoeffs();
131 }
132 
133 
136 (
137  const tmp<scalarField>&
138 ) const
139 {
140  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
141 
142  // For fixedValue nuTilda patches
144 }
145 
146 
149 (
150  const tmp<scalarField>&
151 ) const
152 {
153  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
154 
155  // For zeroGradient nuTilda patches
156  return tmp<Field<scalar>>::New(pos(phip)*(*this));
157 }
158 
159 
161 {
163  writeEntry("value", os);
164  os.writeEntry("solverName", adjointSolverName_);
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
171 (
174 );
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace Foam
179 
180 // ************************************************************************* //
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
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
Foam::adjointFarFieldNuaTildaFvPatchScalarField
Definition: adjointFarFieldNuaTildaFvPatchScalarField.H:53
surfaceFields.H
Foam::surfaceFields.
Foam::adjointBoundaryCondition::boundaryContrPtr_
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Definition: adjointBoundaryCondition.H:73
fvPatchFieldMapper.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::adjointScalarBoundaryCondition
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
Definition: adjointBoundaryConditionsFwd.H:41
adjointFarFieldNuaTildaFvPatchScalarField.H
Foam::Field< vector >
Foam::adjointFarFieldNuaTildaFvPatchScalarField::valueBoundaryCoeffs
virtual tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Definition: adjointFarFieldNuaTildaFvPatchScalarField.C:149
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::adjointFarFieldNuaTildaFvPatchScalarField::adjointFarFieldNuaTildaFvPatchScalarField
adjointFarFieldNuaTildaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: adjointFarFieldNuaTildaFvPatchScalarField.C:45
dict
dictionary dict
Definition: searchingEngine.H:14
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::adjointFarFieldNuaTildaFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: adjointFarFieldNuaTildaFvPatchScalarField.C:160
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::adjointFarFieldNuaTildaFvPatchScalarField::valueInternalCoeffs
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
Definition: adjointFarFieldNuaTildaFvPatchScalarField.C:136
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::adjointFarFieldNuaTildaFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: adjointFarFieldNuaTildaFvPatchScalarField.C:101
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
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::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
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