ATCstandard.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-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 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 #include "ATCstandard.H"
32 #include "wallFvPatch.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(ATCstandard, 0);
43 (
44  ATCModel,
45  ATCstandard,
47 );
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 ATCstandard::ATCstandard
53 (
54  const fvMesh& mesh,
55  const incompressibleVars& primalVars,
56  const incompressibleAdjointVars& adjointVars,
57  const dictionary& dict
58 )
59 :
60  ATCModel(mesh, primalVars, adjointVars, dict),
61  gradU_
62  (
63  IOobject
64  (
65  "gradUATC",
66  mesh_.time().timeName(),
67  mesh_,
68  IOobject::NO_READ,
69  IOobject::NO_WRITE
70  ),
71  mesh_,
73  )
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  addProfiling(ATCstandard, "ATCstandard::addATC");
82  const volVectorField& U = primalVars_.U();
83  const volVectorField& Ua = adjointVars_.UaInst();
85 
86 
87  // Main ATC term
88  ATC_ = gradU_ & Ua;
89 
90  if (extraConvection_ > 0)
91  {
92  // Implicit part added to increase diagonal dominance
94 
95  // Correct rhs due to implicitly augmenting the adjoint convection
96  ATC_ += extraConvection_*(fvc::grad(Ua, "gradUaATC")().T() & U);
97  }
98 
99  //zero ATC on cells next to given patch types
101 
102  // actual ATC term
103  UaEqn += ATC_.internalField();
104 }
105 
106 
108 {
109  const volVectorField& U = primalVars_.U();
110  const volVectorField& Ua = adjointVars_.Ua();
111 
112  // We only need to modify the boundaryField of gradU locally.
113  // If grad(U) is cached then
114  // a. The .ref() call fails since the tmp is initialised from a
115  // const ref
116  // b. we would be changing grad(U) for all other places in the code
117  // that need it
118  // So, always allocate new memory and avoid registering the new field
119  tmp<volTensorField> tgradU(volTensorField::New("gradULocal", fvc::grad(U)));
120  volTensorField::Boundary& gradUbf = tgradU.ref().boundaryFieldRef();
121 
122  // Explicitly correct the boundary gradient to get rid of
123  // the tangential component
124  forAll(mesh_.boundary(), patchI)
125  {
126  const fvPatch& patch = mesh_.boundary()[patchI];
127  if (isA<wallFvPatch>(patch))
128  {
129  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
130  gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
131  }
132  }
133 
134  const volScalarField& mask = getLimiter();
135 
136  return
138  (
139  "ATCFISensitivityTerm" + type(),
140  - (tgradU & Ua)*U*mask
141  );
142 }
143 
144 
146 {
147  const volVectorField& U = primalVars_.U();
149  // Build U to go into the ATC term, based on whether to smooth field or not
150  autoPtr<volVectorField> UForATC(nullptr);
152  {
153  gradU_ = fvc::grad(fvc::reconstruct(phi), "gradUATC");
154  }
155  else
156  {
157  gradU_ = fvc::grad(U, "gradUATC");
158  }
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace Foam
165 
166 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
volScalarField ATClimiter_
Definition: ATCModel.H:91
virtual void addATC(fvVectorMatrix &UaEqn)
Add ATC.
Definition: ATCstandard.C:72
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:44
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Generic dimensioned Type class.
The ATC formualtion resulting by differentiating the Non-conservative form of the momentum equations...
Definition: ATCstandard.H:50
const dimensionSet dimless
Dimensionless.
const incompressibleVars & primalVars_
Definition: ATCModel.H:80
virtual tmp< volTensorField > getFISensitivityTerm() const
Get the FI sensitivity derivatives term coming from the ATC.
Definition: ATCstandard.C:100
Macros for easy insertion into run-time selection tables.
Class including all adjoint fields for incompressible flows.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const incompressibleAdjointVars & adjointVars_
Definition: ATCModel.H:81
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const surfaceScalarField & phi() const
Return const reference to volume flux.
Base class for solution control classes.
#define addProfiling(Name, Descr)
Define profiling trigger with specified name and description string.
dynamicFvMesh & mesh
const fvMesh & mesh_
Definition: ATCModel.H:79
const volVectorField & Ua() const
Return const reference to velocity.
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:56
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
volVectorField ATC_
Definition: ATCModel.H:92
const scalar extraConvection_
Definition: ATCModel.H:85
const volVectorField & U() const
Return const reference to velocity.
static tmp< GeometricField< tensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< tensor >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const bool reconstructGradients_
Definition: ATCModel.H:88
const volVectorField & UaInst() const
Return const reference to velocity.
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:174
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual void updatePrimalBasedQuantities()
Update quantities related with the primal fields.
Definition: ATCstandard.C:138
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127