Go to the documentation of this file.
39 namespace functionObjects
57 if (foundObject<volScalarField>(fieldName_))
63 const volScalarField::Boundary& dpdtBf =
dpdt.boundaryField();
64 const surfaceVectorField::Boundary& SfBf = mesh_.Sf().boundaryField();
68 for (
auto patchi : patchSet_)
70 dfdt.
value() +=
sum(dpdtBf[patchi]*SfBf[patchi]);
85 mesh_.time().timeName(),
93 auto& pDash = tpDash.ref();
98 return store(resultName_, tpDash);
106 Foam::functionObjects::Curle::Curle
120 setResultName(typeName, fieldName_);
136 patchSet_ = mesh_.boundaryMesh().patchSet(
dict.get<
wordRes>(
"patches"));
138 if (patchSet_.empty())
141 <<
"No patches defined"
148 dict.readEntry(
"c0", c0_);
153 const volVectorField::Boundary& Cbf = mesh_.C().boundaryField();
154 const surfaceScalarField::Boundary& magSfBf =
155 mesh_.magSf().boundaryField();
159 for (
auto patchi : patchSet_)
161 x0_.value() +=
sum(Cbf[patchi]*magSfBf[patchi]);
162 sumMagSf +=
sum(magSfBf[patchi]);
168 x0_.value() /= sumMagSf + ROOTVSMALL;
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
static constexpr const zero Zero
Global zero.
const dimensionSet dimVelocity
virtual ~Curle()
Destructor.
Different types of constants.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimArea(sqr(dimLength))
word name(const complex &c)
Return string representation of complex.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &)
Read the Curle data.
constexpr scalar pi(M_PI)
Base class for field expression function objects.
A List of wordRe with additional matching capabilities.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual bool calc()
Calculate the acoustic pressure field and return true if successful.
static const Vector< scalar > zero
Calculate the first temporal derivative.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Graphite solid properties.
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Report a warning using Foam::Warning.