#include <adjointFarFieldPressureFvPatchScalarField.H>


Public Member Functions | |
| TypeName ("adjointFarFieldPressure") | |
| Runtime type information. | |
| adjointFarFieldPressureFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. | |
| adjointFarFieldPressureFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. | |
| adjointFarFieldPressureFvPatchScalarField (const adjointFarFieldPressureFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given adjointOutletPressureFvPatchScalarField. | |
| adjointFarFieldPressureFvPatchScalarField (const adjointFarFieldPressureFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. | |
| virtual tmp< fvPatchField< scalar > > | clone () const |
| Return a clone. | |
| virtual tmp< fvPatchField< scalar > > | clone (const DimensionedField< scalar, volMesh > &iF) const |
| Clone with an internal field reference. | |
| virtual tmp< Field< scalar > > | snGrad () const |
| Return true if this patch field fixes a value. | |
| virtual tmp< Field< scalar > > | valueInternalCoeffs (const tmp< scalarField > &) const |
| Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchField with given weights. | |
| virtual tmp< Field< scalar > > | valueBoundaryCoeffs (const tmp< scalarField > &) const |
| Return the matrix source coefficients corresponding to the. | |
| virtual tmp< Field< scalar > > | gradientInternalCoeffs () const |
| Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this patchField. | |
| virtual tmp< Field< scalar > > | gradientBoundaryCoeffs () const |
| Return the matrix source coefficients corresponding to the evaluation of the gradient of this patchField. | |
| virtual void | updateCoeffs () |
| Update the coefficients associated with the patch field. | |
| virtual void | write (Ostream &) const |
| Write. | |
| virtual void | operator= (const UList< scalar > &) |
| virtual void | operator= (const fvPatchField< scalar > &) |
| virtual void | operator+= (const fvPatchField< scalar > &) |
| virtual void | operator-= (const fvPatchField< scalar > &) |
| virtual void | operator*= (const fvPatchField< scalar > &) |
| virtual void | operator/= (const fvPatchField< scalar > &) |
| virtual void | operator+= (const Field< scalar > &) |
| virtual void | operator-= (const Field< scalar > &) |
| virtual void | operator*= (const Field< scalar > &) |
| virtual void | operator/= (const Field< scalar > &) |
| virtual void | operator= (const scalar) |
| virtual void | operator+= (const scalar) |
| virtual void | operator-= (const scalar) |
| virtual void | operator*= (const scalar) |
| virtual void | operator/= (const scalar) |
| Public Member Functions inherited from adjointBoundaryCondition< scalar > | |
| TypeName ("adjointBoundaryCondition") | |
| Run-time type information. | |
| adjointBoundaryCondition (const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const word &solverName) | |
| Construct from field and base name. | |
| adjointBoundaryCondition (const adjointBoundaryCondition< scalar > &) | |
| Construct as copy. | |
| virtual | ~adjointBoundaryCondition ()=default |
| Destructor. | |
| const word & | objectiveManagerName () const |
| Return objectiveManager name. | |
| const word & | adjointSolverName () const |
| Return adjointSolverName. | |
| const word & | simulationType () const |
| Return the simulationType. | |
| void | setBoundaryContributionPtr () |
| Set the ptr to the correct boundaryAdjointContribution. | |
| boundaryAdjointContribution & | getBoundaryAdjContribution () |
| Get boundaryContribution. | |
| const ATCModel & | getATC () const |
| ATC type might be useful for a number of BCs. Return here. | |
| virtual tmp< Field< typename Foam::outerProduct< Foam::vector, scalar >::type > > | dxdbMult () const |
| Return contribution to sensitivity derivatives. | |
| virtual void | updatePrimalBasedQuantities () |
| Update the primal based quantities related to the adjoint boundary conditions. | |
Additional Inherited Members | |
| Protected Member Functions inherited from adjointBoundaryCondition< scalar > | |
| tmp< Field< typename Foam::outerProduct< Foam::vector, Type2 >::type > > | computePatchGrad (word name) |
| Get gradient of field on a specific boundary. | |
| bool | addATCUaGradUTerm () |
| Whether to add the extra term from the UaGradU formulation. | |
| Protected Attributes inherited from adjointBoundaryCondition< scalar > | |
| const fvPatch & | patch_ |
| Reference to patch. | |
| word | managerName_ |
| objectiveManager name corresponding to field | |
| word | adjointSolverName_ |
| adjointSolver name corresponding to field | |
| word | simulationType_ |
| simulationType corresponding to field. | |
| autoPtr< boundaryAdjointContribution > | boundaryContrPtr_ |
| Engine to manage contributions of the objective functions to the adjoint boundary conditions. | |
| Switch | addATCUaGradUTerm_ |
| Whether to add the extra term from the UaGradU formulation. | |
Definition at line 50 of file adjointFarFieldPressureFvPatchScalarField.H.
| adjointFarFieldPressureFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct from patch and internal field.
Definition at line 32 of file adjointFarFieldPressureFvPatchScalarField.C.
References p.
Referenced by adjointFarFieldPressureFvPatchScalarField(), adjointFarFieldPressureFvPatchScalarField(), and TypeName().

| adjointFarFieldPressureFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const dictionary & | dict ) |
Construct from patch, internal field and dictionary.
Definition at line 58 of file adjointFarFieldPressureFvPatchScalarField.C.
References dict, IOobjectOption::MUST_READ, and p.
| adjointFarFieldPressureFvPatchScalarField | ( | const adjointFarFieldPressureFvPatchScalarField & | ptf, |
| const fvPatch & | p, | ||
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const fvPatchFieldMapper & | mapper ) |
Construct by mapping given adjointOutletPressureFvPatchScalarField.
onto a new patch
Definition at line 44 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointFarFieldPressureFvPatchScalarField(), adjointBoundaryCondition< scalar >::adjointSolverName_, and p.

| adjointFarFieldPressureFvPatchScalarField | ( | const adjointFarFieldPressureFvPatchScalarField & | tppsf, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct as copy setting internal field reference.
Definition at line 73 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointFarFieldPressureFvPatchScalarField().

| TypeName | ( | "adjointFarFieldPressure" | ) |
Runtime type information.
References adjointFarFieldPressureFvPatchScalarField().

|
inlinevirtual |
Return a clone.
Definition at line 109 of file adjointFarFieldPressureFvPatchScalarField.H.
References fvPatchField< Type >::Clone().

|
inlinevirtual |
Clone with an internal field reference.
Definition at line 117 of file adjointFarFieldPressureFvPatchScalarField.H.
References fvPatchField< Type >::Clone().

|
virtual |
Return true if this patch field fixes a value.
Needed to check if a level has to be specified while solving Poissons equations.
Return gradient at boundary
Definition at line 161 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, and Foam::pos().

|
virtual |
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchField with given weights.
Definition at line 176 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, and Foam::neg().

|
virtual |
Return the matrix source coefficients corresponding to the.
evaluation of the value of this patchField with given weights
Definition at line 194 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, and Foam::pos().

|
virtual |
Return the matrix diagonal coefficients corresponding to the evaluation of the gradient of this patchField.
Definition at line 212 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, and Foam::pos().

|
virtual |
Return the matrix source coefficients corresponding to the evaluation of the gradient of this patchField.
Definition at line 228 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, and Foam::pos().

|
virtual |
Update the coefficients associated with the patch field.
Definition at line 87 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::addATCUaGradUTerm(), adjointBoundaryCondition< scalar >::boundaryContrPtr_, delta, Foam::neg(), fvPatchField< Type >::patchInternalField(), Foam::pos(), and tmp< T >::ref().

|
virtual |
Write.
Definition at line 243 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::adjointSolverName_, os(), fvPatchField< Type >::write(), and fvPatchField< Type >::writeValueEntry().

|
virtual |
Definition at line 253 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 265 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::check(), Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 278 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::check(), Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 291 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::check(), Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 304 of file adjointFarFieldPressureFvPatchScalarField.C.
References Foam::abort(), adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::FatalError, FatalErrorInFunction, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 323 of file adjointFarFieldPressureFvPatchScalarField.C.
References Foam::abort(), adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::FatalError, FatalErrorInFunction, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 342 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 354 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
|
virtual |
|
virtual |
Definition at line 390 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 402 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 414 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), and Foam::pos().

|
virtual |
Definition at line 430 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), Foam::pos(), and s().

|
virtual |
Definition at line 442 of file adjointFarFieldPressureFvPatchScalarField.C.
References adjointBoundaryCondition< scalar >::boundaryContrPtr_, Foam::neg(), Field< Type >::operator=(), Foam::pos(), and s().
