Go to the documentation of this file.
42 Foam::timeVaryingMassSorptionFvPatchScalarField::ddtSchemeTypeNames_
45 ddtSchemeType::tsEuler,
46 fv::EulerDdtScheme<scalar>::typeName_()
49 ddtSchemeType::tsCrankNicolson,
50 fv::CrankNicolsonDdtScheme<scalar>::typeName_()
53 ddtSchemeType::tsBackward,
54 fv::backwardDdtScheme<scalar>::typeName_()
68 fixedValueFvPatchScalarField(
p, iF),
83 fixedValueFvPatchScalarField(
p, iF,
dict, false),
88 if (
dict.found(
"value"))
105 const timeVaryingMassSorptionFvPatchScalarField& ptf,
107 const DimensionedField<scalar, volMesh>& iF,
108 const fvPatchFieldMapper& mapper
111 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
124 fixedValueFvPatchScalarField(ptf),
138 fixedValueFvPatchScalarField(ptf, iF),
152 fixedValueFvPatchScalarField::autoMap(m);
162 fixedValueFvPatchScalarField::rmap(ptf, addr);
170 auto& source = tsource.ref();
175 source = -kabs_*w*
max(patchInternalField() -
cp, scalar(0));
177 source += kdes_*
max(
cp - patchInternalField(), scalar(0));
190 const label patchi =
patch().index();
192 const scalar dt = db().time().deltaTValue();
199 const word ddtSchemeName(
fld.mesh().ddtScheme(
fld.name()));
200 const ddtSchemeType& ddtScheme = ddtSchemeTypeNames_[ddtSchemeName];
205 tmp<scalarField> dfldp =
208 *
max(patchInternalField() -
cp, scalar(0))
213 *
max(
cp - patchInternalField(), scalar(0))
219 case tsCrankNicolson:
221 operator==(fld0.boundaryField()[patchi] + dfldp);
227 const scalar dt0 = db().time().deltaT0Value();
229 const scalar
c = scalar(1) + dt/(dt + dt0);
230 const scalar c00 = dt*dt/(dt0*(dt + dt0));
231 const scalar c0 =
c + c00;
236 c0*fld0.boundaryField()[patchi]
237 - c00*fld0.oldTime().boundaryField()[patchi]
247 << ddtSchemeName <<
nl
248 <<
" on patch " << this->
patch().name()
249 <<
" of field " << this->internalField().name()
250 <<
" in file " << this->internalField().objectPath()
255 fixedValueFvPatchScalarField::updateCoeffs();
267 writeEntry(
"value",
os);
278 timeVaryingMassSorptionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
virtual void operator=(const UList< Type > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void write(Ostream &) const
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
timeVaryingMassSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual void autoMap(const fvPatchFieldMapper &)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual void rmap(const fvPatchScalarField &, const labelList &)
static tmp< T > New(Args &&... args)
Ostream & writeEntry(const keyType &key, const T &value)
const dimensionedScalar c
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition provides a first order fixed-value condition for a given scalar field to mode...
Generic GeometricField class.
tmp< scalarField > source() const
A min/max value pair with additional methods. In addition to conveniently storing values,...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...