Go to the documentation of this file.
39 namespace incompressible
79 for (
const label patchI : sensitivityPatchIDs_)
91 Info<<
"Mesh Movement Propagation(direction, CP), ("
92 << idir <<
", " << iCP <<
"), Iteration : "<< iter <<
endl;
100 scalar residual =
mag(mEqn.solve().initialResidual());
109 Info<<
"\n***Reached dxdb convergence limit, iteration " << iter
121 sensitivityBezierFI::sensitivityBezierFI
140 flowSens_(3*Bezier_.nBezier(),
Zero),
141 dSdbSens_(3*Bezier_.nBezier(),
Zero),
142 dndbSens_(3*Bezier_.nBezier(),
Zero),
143 dxdbDirectSens_(3*Bezier_.nBezier(),
Zero),
144 dVdbSens_(3*Bezier_.nBezier(),
Zero),
145 distanceSens_(3*Bezier_.nBezier(),
Zero),
146 optionsSens_(3*Bezier_.nBezier(),
Zero),
147 bcSens_(3*Bezier_.nBezier(),
Zero),
149 derivativesFolder_(
"optimisation"/
type() +
"Derivatives"),
151 meshMovementIters_(-1),
152 meshMovementResidualLimit_(1.
e-7),
183 distanceSensPtr.reset
185 createZeroFieldPtr<tensor>
196 const label nDVs = 3*nBezier;
197 for (label iDV = 0; iDV < nDVs; iDV++)
199 label iCP = iDV%nBezier;
200 label idir = iDV/nBezier;
226 for (
const label patchI : sensitivityPatchIDs_)
232 tmp<vectorField> tdSdb =
235 tmp<vectorField> tdndb =
246 tmp<vectorField> tdxdbFace =
275 distanceSensPtr().primitiveField()
276 && gradDxDb.primitiveField()
318 Info<<
"Writing control point sensitivities to file" <<
endl;
329 <<
setw(widthDV) <<
"#dv" <<
" "
330 <<
setw(width) <<
"total" <<
" "
331 <<
setw(width) <<
"flow" <<
" "
332 <<
setw(width) <<
"dSdb" <<
" "
333 <<
setw(width) <<
"dndb" <<
" "
334 <<
setw(width) <<
"dxdbDirect" <<
" "
335 <<
setw(width) <<
"dVdb" <<
" "
336 <<
setw(width) <<
"distance" <<
" "
337 <<
setw(width) <<
"options" <<
" "
342 label lastActive(-1);
344 for (label iDV = 0; iDV < nDVs; iDV++)
346 const label iCP(iDV%nBezier);
347 const label idir(iDV/nBezier);
348 if (!confineMovement[idir][iCP])
350 if (iDV!=lastActive + 1)
356 <<
setw(widthDV) << iDV <<
" "
vectorField optionsDxDbMult_
virtual void clearSensitivities()
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
class for managing incompressible objective functions.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
const word & solverName() const
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
autoPtr< boundaryVectorField > bcDxDbMult_
const boolList & confineXmovement() const
A class for handling words, derived from Foam::string.
virtual void clearSensitivities()
defineTypeNameAndDebug(adjointEikonalSolver, 0)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
const Internal::FieldType & primitiveField() const
bool read(const char *buf, int32_t &val)
volTensorField gradDxDbMult_
static word timeName(const scalar t, const int precision=precision_)
scalarField dxdbDirectSens_
scalarField distanceSens_
static bool master(const label communicator=worldComm)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Ostream & endl(Ostream &os)
Type gSum(const FieldField< Field, Type > &f)
const boolList & confineZmovement() const
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
virtual void assembleSensitivities()
Class including all adjoint fields for incompressible flows.
Field< vector > vectorField
Specialisation of Field<T> for vector.
autoPtr< adjointEikonalSolver > eikonalSolver_
Base class for creating a set of variables.
SolverPerformance< Type > solve(const dictionary &)
Ostream & printExecutionTime(OSstream &os) const
tmp< volVectorField > solveMeshMovementEqn(const label iCP, const label idir)
autoPtr< boundaryVectorField > dSfdbMult_
addToRunTimeSelectionTable(adjointSensitivity, sensitivityBezier, dictionary)
Abstract base class for adjoint-based sensitivities in incompressible flows.
Generic templated field type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Istream and Ostream manipulators taking arguments.
label max(const labelHashSet &set, label maxValue=labelMin)
tmp< tensorField > dxdbFace(const label patchI, const label cpI, bool useChainRule=true) const
static tmp< volVectorField > autoCreateMeshMovementField(const fvMesh &mesh, const word &name, const dimensionSet &dims)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
autoPtr< boundaryVectorField > dxdbDirectMult_
Macros for easy insertion into run-time selection tables.
incompressibleAdjointVars & adjointVars_
Mesh data needed to do the Finite Volume discretisation.
Omanip< int > setw(const int i)
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
scalar meshMovementResidualLimit_
Output to file stream, using an OSstream.
GeometricField< vector, fvPatchField, volMesh > volVectorField
void reset(autoPtr< T > &&other) noexcept
static unsigned int defaultPrecision() noexcept
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Calculation of adjoint based sensitivities for Bezier control points using the FI appoach.
Base class for Field Integral-based sensitivity derivatives.
Internal & ref(const bool updateAccessTime=true)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
const dimensionedScalar e
A special matrix type and solver, designed for finite volume solutions of scalar equations....
word name(const expressions::valueTypeCode typeCode)
const Time & time() const
const boolList & confineYmovement() const
const boolListList & confineMovement() const
virtual void write(const word &baseName=word::null)
autoPtr< boundaryVectorField > dnfdbMult_
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Generic GeometricField class.
bool mkDir(const fileName &pathName, mode_t mode=0777)
Type gMax(const FieldField< Field, Type > &f)
Base class for solution control classes.
fileName derivativesFolder_
const DimensionedField< scalar, volMesh > & V() const
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool returnDimensionedNormalSens=true) const