solveBulkSurfactant.H
Go to the documentation of this file.
1 {
2  volScalarField& C = bulkSurfactantConcentration();
3 
4  const dimensionedScalar& D = surfactant().bulkDiffusion();
5 
6  scalar ka = surfactant().adsorptionCoeff().value();
7  scalar kb = surfactant().desorptionCoeff().value();
8  scalar CsInf = surfactant().saturatedConc().value();
9 
10  const scalarField& Cs =
11  surfactantConcentration().internalField();
12 
13  // const scalarField& Cfs = C.boundaryField()[fsPatchIndex()];
14 
15  if
16  (
17  C.boundaryField()[fsPatchIndex()].type()
18  == fixedGradientFvPatchScalarField::typeName
19  )
20  {
21  fixedGradientFvPatchScalarField& fsC =
22  refCast<fixedGradientFvPatchScalarField>
23  (
24  C.boundaryFieldRef()[fsPatchIndex()]
25  );
26 
27  fsC.gradient() = (ka*kb*Cs - ka*(CsInf-Cs)*fsC)/D.value();
28  }
29  else
30  {
32  << "Bulk concentration boundary condition "
33  << "at the free-surface patch is not "
34  << fixedGradientFvPatchScalarField::typeName
35  << exit(FatalError);
36  }
37 
39  (
40  fvm::ddt(C)
41  + fvm::div(phi(), C, "div(phi,C)")
42  - fvm::laplacian(D, C, "laplacian(D,C)")
43  );
44 
45  CEqn.solve();
46 }
CsInf
scalar CsInf
Definition: solveBulkSurfactant.H:8
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
CEqn
fvScalarMatrix CEqn(fvm::ddt(C)+fvm::div(phi(), C, "div(phi,C)") - fvm::laplacian(D, C, "laplacian(D,C)"))
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::FatalError
error FatalError
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Cs
const scalarField & Cs
Definition: solveBulkSurfactant.H:10
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
kb
scalar kb
Definition: solveBulkSurfactant.H:7
ka
scalar ka
Definition: solveBulkSurfactant.H:6