Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
30#include "ATCstandard.H"
32#include "wallFvPatch.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
43(
45 ATCstandard,
47);
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52ATCstandard::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 (
64 (
65 "gradUATC",
67 mesh_,
70 ),
71 mesh_,
73 )
74{}
75
76
77// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78
80{
81 addProfiling(ATCstandard, "ATCstandard::addATC");
82 const volVectorField& U = primalVars_.U();
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
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition ATCModel.H:59
const incompressibleAdjointVars & adjointVars_
Definition ATCModel.H:81
const fvMesh & mesh_
Definition ATCModel.H:79
const scalar extraConvection_
Definition ATCModel.H:85
volVectorField ATC_
Definition ATCModel.H:92
const bool reconstructGradients_
Definition ATCModel.H:88
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition ATCModel.C:174
void smoothATC()
Limit ATC field using ATClimiter_.
Definition ATCModel.C:44
const incompressibleVars & primalVars_
Definition ATCModel.H:80
volScalarField ATClimiter_
Definition ATCModel.H:91
The ATC formualtion resulting by differentiating the Non-conservative form of the momentum equations.
Definition ATCstandard.H:53
virtual void updatePrimalBasedQuantities()
Update quantities related with the primal fields.
virtual tmp< volTensorField > getFISensitivityTerm() const
Get the FI sensitivity derivatives term coming from the ATC.
virtual void addATC(fvVectorMatrix &UaEqn)
Add ATC.
Definition ATCstandard.C:72
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())
GeometricBoundaryField< tensor, fvPatchField, volMesh > Boundary
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const volVectorField & Ua() const
Return const reference to velocity.
const volVectorField & UaInst() const
Return const reference to velocity.
Class including all adjoint fields for incompressible flows.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
const surfaceScalarField & phi() const
Return const reference to volume flux.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & T
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299