Go to the documentation of this file.
57 bool nbrModelFound =
false;
65 refCast<const interRegionHeatTransferModel>(
fvOptions[i])
76 <<
" in region " << nbrMesh.
name() <<
nl
91 if (mesh_.time().timeIndex() != timeIndex_)
94 timeIndex_ = mesh_.time().timeIndex();
110 const word& modelType,
122 nbrModelName_(coeffs_.lookup(
"nbrModelName")),
131 mesh.time().timeName(),
143 zeroGradientFvPatchScalarField::typeName
145 semiImplicit_(
false),
146 TName_(coeffs_.lookupOrDefault<
word>(
"TName",
"T")),
147 TNbrName_(coeffs_.lookupOrDefault<
word>(
"TNbrName",
"T"))
151 coeffs_.lookup(
"fieldNames") >> fieldNames_;
152 applied_.setSize(fieldNames_.size(),
false);
154 coeffs_.lookup(
"semiImplicit") >> semiImplicit_;
188 mesh_.time().timeName(),
208 Info<<
"Volumetric integral of htc: "
212 if (mesh_.time().outputTime())
230 eqn += htc_*(Tmapped -
T) + htcByCpv*
he -
fvm::Sp(htcByCpv,
he);
237 Info<<
"Energy exchange from region " << nbrMesh.
name()
238 <<
" To " << mesh_.name() <<
" : " << energy.
value()
245 <<
" on mesh " << mesh_.name()
246 <<
" could not find object basicThermo."
247 <<
" The available objects are: "
259 eqn += htc_*(Tmapped -
T);
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
word nbrRegionName_
Name of the neighbour region to map.
void correct()
Correct to calculate the inter-region heat transfer coefficient.
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionSet dimEnergy
const word & name() const
Return const access to the source name.
static const word dictName
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffis...
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldI)
Source term to energy equation.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Abstract base-class for fluid and solid thermodynamic properties.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
interRegionHeatTransferModel * nbrModel_
Pointer to neighbour interRegionHeatTransferModel.
const fvMesh & mesh_
Reference to the mesh database.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Base class for inter-region exchange.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word nbrModelName_
Name of the model in the neighbour mesh.
InternalField & internalField()
Return internal field.
interRegionHeatTransferModel(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
bool firstIter_
First iteration.
Calculate the matrix for implicit and explicit sources.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void setNbrModel()
Set the neighbour interRegionHeatTransferModel.
defineTypeNameAndDebug(option, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual ~interRegionHeatTransferModel()
Destructor.
Volume integrate volField creating a volField.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Time & time() const
Return the top-level database.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
const dimensionSet dimVolume(pow3(dimLength))
Generic GeometricField class.
List of finite volume options.
word name(const complex &)
Return a string representation of a complex.
const word & name() const
Return reference to name.