/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . \*---------------------------------------------------------------------------*/ #include "multiphaseSystem.H" #include "alphaContactAngleFvPatchScalarField.H" #include "Time.H" #include "subCycle.H" #include "MULES.H" #include "fvcSnGrad.H" #include "fvcFlux.H" // * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * // const Foam::scalar Foam::multiphaseSystem::convertToRad = Foam::constant::mathematical::pi/180.0; // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // void Foam::multiphaseSystem::calcAlphas() { scalar level = 0.0; alphas_ == 0.0; forAllIter(PtrDictionary, phases_, iter) { alphas_ += level*iter(); level += 1.0; } alphas_.correctBoundaryConditions(); } void Foam::multiphaseSystem::solveAlphas() { word alphaScheme("div(phi,alpha)"); word alpharScheme("div(phirb,alpha)"); surfaceScalarField phic(mag(phi_/mesh_.magSf())); PtrList phiAlphaCorrs(phases_.size()); int phasei = 0; forAllIter(PtrDictionary, phases_, iter) { phaseModel& phase1 = iter(); phase1.phiAlpha() = dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0); phiAlphaCorrs.set ( phasei, new surfaceScalarField ( fvc::flux ( phi_, phase1, alphaScheme ) ) ); surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei]; forAllIter(PtrDictionary, phases_, iter2) { phaseModel& phase2 = iter2(); if (&phase2 == &phase1) continue; surfaceScalarField phir ( (phase1.phi() - phase2.phi()) + min(cAlpha(phase1, phase2)*phic, max(phic)) *nHatf(phase1, phase2) ); phiAlphaCorr += fvc::flux ( -fvc::flux(-phir, phase2, alpharScheme), phase1, alpharScheme ); } MULES::limit ( geometricOneField(), phase1, phi_, phiAlphaCorr, zeroField(), zeroField(), 1, 0, 3, true ); phasei++; } MULES::limitSum(phiAlphaCorrs); volScalarField sumAlpha ( IOobject ( "sumAlpha", mesh_.time().timeName(), mesh_ ), mesh_, dimensionedScalar("sumAlpha", dimless, 0) ); phasei = 0; forAllIter(PtrDictionary, phases_, iter) { phaseModel& phase1 = iter(); surfaceScalarField& phiAlpha = phiAlphaCorrs[phasei]; phiAlpha += upwind(mesh_, phi_).flux(phase1); MULES::explicitSolve ( geometricOneField(), phase1, phiAlpha, zeroField(), zeroField() ); phase1.phiAlpha() += phiAlpha; Info<< phase1.name() << " volume fraction, min, max = " << phase1.weightedAverage(mesh_.V()).value() << ' ' << min(phase1).value() << ' ' << max(phase1).value() << endl; sumAlpha += phase1; phasei++; } Info<< "Phase-sum volume fraction, min, max = " << sumAlpha.weightedAverage(mesh_.V()).value() << ' ' << min(sumAlpha).value() << ' ' << max(sumAlpha).value() << endl; calcAlphas(); } Foam::dimensionedScalar Foam::multiphaseSystem::sigma ( const phaseModel& phase1, const phaseModel& phase2 ) const { scalarCoeffTable::const_iterator sigma ( sigmas_.find(interfacePair(phase1, phase2)) ); if (sigma == sigmas_.end()) { FatalErrorIn ( "multiphaseSystem::sigma(const phaseModel& phase1," "const phaseModel& phase2) const" ) << "Cannot find interface " << interfacePair(phase1, phase2) << " in list of sigma values" << exit(FatalError); } return dimensionedScalar("sigma", dimSigma_, sigma()); } Foam::scalar Foam::multiphaseSystem::cAlpha ( const phaseModel& phase1, const phaseModel& phase2 ) const { scalarCoeffTable::const_iterator cAlpha ( cAlphas_.find(interfacePair(phase1, phase2)) ); if (cAlpha == cAlphas_.end()) { FatalErrorIn ( "multiphaseSystem::cAlpha" "(const phaseModel& phase1, const phaseModel& phase2) const" ) << "Cannot find interface " << interfacePair(phase1, phase2) << " in list of cAlpha values" << exit(FatalError); } return cAlpha(); } Foam::dimensionedScalar Foam::multiphaseSystem::Cvm ( const phaseModel& phase1, const phaseModel& phase2 ) const { scalarCoeffTable::const_iterator Cvm ( Cvms_.find(interfacePair(phase1, phase2)) ); if (Cvm != Cvms_.end()) { return Cvm()*phase2.rho(); } Cvm = Cvms_.find(interfacePair(phase2, phase1)); if (Cvm != Cvms_.end()) { return Cvm()*phase1.rho(); } FatalErrorIn ( "multiphaseSystem::sigma" "(const phaseModel& phase1, const phaseModel& phase2) const" ) << "Cannot find interface " << interfacePair(phase1, phase2) << " in list of sigma values" << exit(FatalError); return Cvm()*phase2.rho(); } Foam::tmp Foam::multiphaseSystem::nHatfv ( const volScalarField& alpha1, const volScalarField& alpha2 ) const { /* // Cell gradient of alpha volVectorField gradAlpha = alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2); // Interpolated face-gradient of alpha surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha); */ surfaceVectorField gradAlphaf ( fvc::interpolate(alpha2)*fvc::interpolate(fvc::grad(alpha1)) - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2)) ); // Face unit interface normal return gradAlphaf/(mag(gradAlphaf) + deltaN_); } Foam::tmp Foam::multiphaseSystem::nHatf ( const volScalarField& alpha1, const volScalarField& alpha2 ) const { // Face unit interface normal flux return nHatfv(alpha1, alpha2) & mesh_.Sf(); } // Correction for the boundary condition on the unit normal nHat on // walls to produce the correct contact angle. // The dynamic contact angle is calculated from the component of the // velocity on the direction of the interface, parallel to the wall. void Foam::multiphaseSystem::correctContactAngle ( const phaseModel& phase1, const phaseModel& phase2, surfaceVectorField::GeometricBoundaryField& nHatb ) const { const volScalarField::GeometricBoundaryField& gbf = phase1.boundaryField(); const fvBoundaryMesh& boundary = mesh_.boundary(); forAll(boundary, patchi) { if (isA(gbf[patchi])) { const alphaContactAngleFvPatchScalarField& acap = refCast(gbf[patchi]); vectorField& nHatPatch = nHatb[patchi]; vectorField AfHatPatch ( mesh_.Sf().boundaryField()[patchi] /mesh_.magSf().boundaryField()[patchi] ); alphaContactAngleFvPatchScalarField::thetaPropsTable:: const_iterator tp = acap.thetaProps().find(interfacePair(phase1, phase2)); if (tp == acap.thetaProps().end()) { FatalErrorIn ( "multiphaseSystem::correctContactAngle" "(const phaseModel& phase1, const phaseModel& phase2, " "fvPatchVectorFieldField& nHatb) const" ) << "Cannot find interface " << interfacePair(phase1, phase2) << "\n in table of theta properties for patch " << acap.patch().name() << exit(FatalError); } bool matched = (tp.key().first() == phase1.name()); scalar theta0 = convertToRad*tp().theta0(matched); scalarField theta(boundary[patchi].size(), theta0); scalar uTheta = tp().uTheta(); // Calculate the dynamic contact angle if required if (uTheta > SMALL) { scalar thetaA = convertToRad*tp().thetaA(matched); scalar thetaR = convertToRad*tp().thetaR(matched); // Calculated the component of the velocity parallel to the wall vectorField Uwall ( phase1.U().boundaryField()[patchi].patchInternalField() - phase1.U().boundaryField()[patchi] ); Uwall -= (AfHatPatch & Uwall)*AfHatPatch; // Find the direction of the interface parallel to the wall vectorField nWall ( nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch ); // Normalise nWall nWall /= (mag(nWall) + SMALL); // Calculate Uwall resolved normal to the interface parallel to // the interface scalarField uwall(nWall & Uwall); theta += (thetaA - thetaR)*tanh(uwall/uTheta); } // Reset nHatPatch to correspond to the contact angle scalarField a12(nHatPatch & AfHatPatch); scalarField b1(cos(theta)); scalarField b2(nHatPatch.size()); forAll(b2, facei) { b2[facei] = cos(acos(a12[facei]) - theta[facei]); } scalarField det(1.0 - a12*a12); scalarField a((b1 - a12*b2)/det); scalarField b((b2 - a12*b1)/det); nHatPatch = a*AfHatPatch + b*nHatPatch; nHatPatch /= (mag(nHatPatch) + deltaN_.value()); } } } Foam::tmp Foam::multiphaseSystem::K ( const phaseModel& phase1, const phaseModel& phase2 ) const { tmp tnHatfv = nHatfv(phase1, phase2); correctContactAngle(phase1, phase2, tnHatfv().boundaryField()); // Simple expression for curvature return -fvc::div(tnHatfv & mesh_.Sf()); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::multiphaseSystem::multiphaseSystem ( const fvMesh& mesh, const surfaceScalarField& phi ) : IOdictionary ( IOobject ( "transportProperties", mesh.time().constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE ) ), phases_(lookup("phases"), phaseModel::iNew(mesh)), mesh_(mesh), phi_(phi), alphas_ ( IOobject ( "alphas", mesh_.time().timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), mesh_, dimensionedScalar("alphas", dimless, 0.0), zeroGradientFvPatchScalarField::typeName ), sigmas_(lookup("sigmas")), dimSigma_(1, 0, -2, 0, 0), cAlphas_(lookup("interfaceCompression")), Cvms_(lookup("virtualMass")), deltaN_ ( "deltaN", 1e-8/pow(average(mesh_.V()), 1.0/3.0) ) { calcAlphas(); alphas_.write(); forAllIter(PtrDictionary, phases_, iter) { phaseModelTable_.add(iter()); } interfaceDictTable dragModelsDict(lookup("drag")); forAllConstIter(interfaceDictTable, dragModelsDict, iter) { dragModels_.insert ( iter.key(), dragModel::New ( iter(), *phaseModelTable_.find(iter.key().first())(), *phaseModelTable_.find(iter.key().second())() ).ptr() ); } } // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // Foam::tmp Foam::multiphaseSystem::rho() const { PtrDictionary::const_iterator iter = phases_.begin(); tmp trho = iter()*iter().rho(); for (++iter; iter != phases_.end(); ++iter) { trho() += iter()*iter().rho(); } return trho; } Foam::tmp Foam::multiphaseSystem::Cvm ( const phaseModel& phase ) const { tmp tCvm ( new volScalarField ( IOobject ( "Cvm", mesh_.time().timeName(), mesh_ ), mesh_, dimensionedScalar ( "Cvm", dimensionSet(1, -3, 0, 0, 0), 0 ) ) ); forAllConstIter(PtrDictionary, phases_, iter) { const phaseModel& phase2 = iter(); if (&phase2 != &phase) { tCvm() += Cvm(phase, phase2)*phase2; } } return tCvm; } Foam::tmp Foam::multiphaseSystem::Svm ( const phaseModel& phase ) const { tmp tSvm ( new volVectorField ( IOobject ( "Svm", mesh_.time().timeName(), mesh_ ), mesh_, dimensionedVector ( "Svm", dimensionSet(1, -2, -2, 0, 0), vector::zero ) ) ); forAllConstIter(PtrDictionary, phases_, iter) { const phaseModel& phase2 = iter(); if (&phase2 != &phase) { tSvm() += Cvm(phase, phase2)*phase2*phase2.DDtU(); } } return tSvm; } Foam::autoPtr Foam::multiphaseSystem::dragCoeffs() const { autoPtr dragCoeffsPtr(new dragCoeffFields); forAllConstIter(dragModelTable, dragModels_, iter) { const dragModel& dm = *iter(); dragCoeffsPtr().insert ( iter.key(), ( dm.phase1()*dm.phase2() *dm.K(mag(dm.phase1().U() - dm.phase2().U())) + dm.residualDrag() ).ptr() ); } return dragCoeffsPtr; } Foam::tmp Foam::multiphaseSystem::dragCoeff ( const phaseModel& phase, const dragCoeffFields& dragCoeffs ) const { tmp tdragCoeff ( new volScalarField ( IOobject ( "dragCoeff", mesh_.time().timeName(), mesh_ ), mesh_, dimensionedScalar ( "dragCoeff", dimensionSet(1, -3, -1, 0, 0), 0 ) ) ); dragModelTable::const_iterator dmIter = dragModels_.begin(); dragCoeffFields::const_iterator dcIter = dragCoeffs.begin(); for ( ; dmIter != dragModels_.end() && dcIter != dragCoeffs.end(); ++dmIter, ++dcIter ) { if ( &phase == &dmIter()->phase1() || &phase == &dmIter()->phase2() ) { tdragCoeff() += *dcIter(); } } return tdragCoeff; } Foam::tmp Foam::multiphaseSystem::surfaceTensionForce() const { tmp tstf ( new surfaceScalarField ( IOobject ( "surfaceTensionForce", mesh_.time().timeName(), mesh_ ), mesh_, dimensionedScalar ( "surfaceTensionForce", dimensionSet(1, -2, -2, 0, 0), 0.0 ) ) ); surfaceScalarField& stf = tstf(); forAllConstIter(PtrDictionary, phases_, iter1) { const phaseModel& phase1 = iter1(); PtrDictionary::const_iterator iter2 = iter1; ++iter2; for (; iter2 != phases_.end(); ++iter2) { const phaseModel& phase2 = iter2(); stf += sigma(phase1, phase2) *fvc::interpolate(K(phase1, phase2))* ( fvc::interpolate(phase2)*fvc::snGrad(phase1) - fvc::interpolate(phase1)*fvc::snGrad(phase2) ); } } return tstf; } Foam::tmp Foam::multiphaseSystem::nearInterface() const { tmp tnearInt ( new volScalarField ( IOobject ( "nearInterface", mesh_.time().timeName(), mesh_ ), mesh_, dimensionedScalar("nearInterface", dimless, 0.0) ) ); forAllConstIter(PtrDictionary, phases_, iter) { tnearInt() = max(tnearInt(), pos(iter() - 0.01)*pos(0.99 - iter())); } return tnearInt; } void Foam::multiphaseSystem::solve() { forAllIter(PtrDictionary, phases_, iter) { iter().correct(); } const Time& runTime = mesh_.time(); const dictionary& pimpleDict = mesh_.solutionDict().subDict("PIMPLE"); label nAlphaSubCycles(readLabel(pimpleDict.lookup("nAlphaSubCycles"))); volScalarField& alpha = phases_.first(); if (nAlphaSubCycles > 1) { dimensionedScalar totalDeltaT = runTime.deltaT(); for ( subCycle alphaSubCycle(alpha, nAlphaSubCycles); !(++alphaSubCycle).end(); ) { solveAlphas(); } } else { solveAlphas(); } } bool Foam::multiphaseSystem::read() { if (regIOobject::read()) { bool readOK = true; PtrList phaseData(lookup("phases")); label phasei = 0; forAllIter(PtrDictionary, phases_, iter) { readOK &= iter().read(phaseData[phasei++].dict()); } lookup("sigmas") >> sigmas_; lookup("interfaceCompression") >> cAlphas_; lookup("virtualMass") >> Cvms_; return readOK; } else { return false; } } // ************************************************************************* //