|
[Sponsors] |
Adding the Energy Equation to interFoam (OF 2.4.0) |
|
LinkBack | Thread Tools | Search this Thread | Display Modes |
November 13, 2015, 06:35 |
Adding the Energy Equation to interFoam (OF 2.4.0)
|
#1 |
New Member
Perry Uhlig
Join Date: Oct 2015
Location: Germany
Posts: 13
Rep Power: 11 |
Dear Foamers,
first of all, this is my first thread, me being new to OF and the forum as well. So I am thankful for any advice regarding anything related to CFD in general and OF and this forum in particular. Adding the Energy Equation to interFaom makes it necessary to modify the used library (A) and secondly modify the solver itself (B). I am having problems compiling said application, hence I have written this thread. I am going to describe my steps in detail in this thread. The procedure is inspired by both the documents from Damiano Natali ("interTempFoam") and Qingming Liu ("myinterFoamDiabatic"). The links are given beneath. http://www.wolfdynamics.com/images/c...erTempFoam.pdf http://www.tfd.chalmers.se/~hani/kur...gLIU-final.pdf A) Modifcation of the library 1) Copy the original folder transportModels to the user's directory cd $WM_PROJECT_USER_DIR/src/ cp -rp $FOAM_SRC/transportModels . 2) Rename the folders and files which are being modified (namely the libraries incompressibleTwoPhaseMixture and immiscibleIncompressibleTwoPhaseMixture). cd transportModels/incompressible mv incompressibleTwoPhaseMixture myIncompressibleTwoPhaseMixture cd myIncompressibleTwoPhaseMixture mv incompressibleTwoPhaseMixture.C myIncompressibleTwoPhaseMixture.C mv incompressibleTwoPhaseMixture.H myIncompressibleTwoPhaseMixture.H cd ../.. mv immiscibleIncompressibleTwoPhaseMixture myImmiscibleIncompressibleTwoPhaseMixture cd myImmiscibleIncompressibleTwoPhaseMixture mv immiscibleIncompressibleTwoPhaseMixture.C myImmiscibleIncompressibleTwoPhaseMixture.C mv immiscibleIncompressibleTwoPhaseMixture.H myImmiscibleIncompressibleTwoPhaseMixture.H 3) Change the respective file files corresponding to the the naming and path. in the folder incompressible/Make for myIncompressibleTwoPhaseMixture: Code:
myIncompressibleTwoPhaseMixture/myIncompressibleTwoPhaseMixture.C LIB = $(FOAM_USER_LIBBIN)/libmyIncompressibleTransportModels Code:
myImmiscibleIncompressibleTwoPhaseMixture.C LIB = $(FOAM_USER_LIBBIN)/libmyImmiscibleIncompressibleTwoPhaseMixture in myIncompressibleTwoPhaseMixture.C: Code:
#include "myIncompressibleTwoPhaseMixture.H" Code:
#include "myImmiscibleIncompressibleTwoPhaseMixture.H" Code:
#include "myIncompressibleTwoPhaseMixture.H" In the file myImmiscibleIncompressibleTwoPhaseMixture/Make/options make the adaption as highlighted, so the connection to the new library is made. Code:
LIB_LIBS = \ -ltwoPhaseMixture \ -L$(FOAM_USER_LIBBIN) \ -lmyIncompressibleTransportModels \ -linterfaceProperties \ -ltwoPhaseProperties \ -lfiniteVolume Add the new parameters rho and Pr in myIncompressibleTwoPhaseMixture.H. Code:
dimensionedScalar rho1_; dimensionedScalar rho2_; //-----------modified_start------------// dimensionedScalar cp1_; dimensionedScalar cp2_; dimensionedScalar Pr1_; dimensionedScalar Pr2_; //-----------modified_end------------// Code:
const dimensionedScalar& rho2() const { return rho2_; }; //-----------modified_start------------// //- Return const-access to phase1 density const dimensionedScalar& cp1() const { return cp1_; } //- Return const-access to phase2 density const dimensionedScalar& cp2() const { return cp2_; }; //- Return const-access to phase1 density const dimensionedScalar& Pr1() const { return Pr1_; } //- Return const-access to phase2 density const dimensionedScalar& Pr2() const { return Pr2_; }; //-----------modified_end------------// Code:
tmp<surfaceScalarField> nuf() const; //-----------modified_start------------// //- Return the face-interpolated conductivity tmp<surfaceScalarField> kappaf() const; //-----------modified_end------------// In myIncompressibleTwoPhaseMixture.C add the operations for the new parameters. Code:
rho2_("rho", dimDensity, nuModel2_->viscosityProperties().lookup("rho")), //-----------modified_start------------// cp1_("cp", dimensionSet(0, 2, -2, -1, 0, 0, 0), nuModel1_->viscosityProperties().lookup("cp")), cp2_("cp", dimensionSet(0, 2, -2, -1, 0, 0, 0), nuModel2_->viscosityProperties().lookup("cp")), Pr1_("Pr", dimensionSet(0, 0, 0, 0, 0, 0, 0), nuModel1_->viscosityProperties().lookup("Pr")), Pr2_("Pr", dimensionSet(0, 0, 0, 0, 0, 0, 0), nuModel2_->viscosityProperties().lookup("Pr")), //-----------modified_end------------// Code:
nuModel2_->viscosityProperties().lookup("rho") >> rho2_; //-----------modified_start--------------------// nuModel1_->viscosityProperties().lookup("cp") >> cp1_; nuModel2_->viscosityProperties().lookup("cp") >> cp2_; nuModel1_->viscosityProperties().lookup("Pr") >> Pr1_; nuModel2_->viscosityProperties().lookup("Pr") >> Pr2_; //-----------modified_end----------------------// Code:
Foam::tmp<Foam::surfaceScalarField> Foam::incompressibleTwoPhaseMixture::nuf() const { const surfaceScalarField alpha1f ( min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1)) ); return tmp<surfaceScalarField> ( new surfaceScalarField ( "nuf", ( alpha1f*rho1_*fvc::interpolate(nuModel1_->nu()) + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu()) )/(alpha1f*rho1_ + (scalar(1) - alpha1f)*rho2_) ) ); } //-----------modified_start------------// Foam::tmp<Foam::surfaceScalarField> Foam::incompressibleTwoPhaseMixture::kappaf() const { const surfaceScalarField alpha1f ( min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1)) ); return tmp<surfaceScalarField> ( new surfaceScalarField ( "kappaf", ( alpha1f*rho1_*cp1_*(1/Pr1_) *fvc::interpolate(nuModel1_->nu()) + (scalar(1) - alpha1f)*rho2_*cp2_ *(1/Pr2_)*fvc::interpolate(nuModel2_->nu()) ) ) ); } //-----------modified_end------------// Compile both libraries from the respective folder. Start with myIncompressibleTwoPhaseMixture, because myImmiscibleIncompressibleTwoPhaseMixture depends on it already being compiled beforehand. cd incompressible wclean wmake libso cd .. cd myImmiscibleIncompressibleTwoPhaseMixture wclean wmake libso Compiling works at this stage, which is good for the moment but doesn't necessarily mean to me that I did everything right. B) Modification of the solver interFoam to interTempFoam 1) Enter the destination of the new solver and copy the base solver. cd $WM_PROJECT_USER_DIR/applications/solvers/multiphase cp -rp $FOAM_APP/solvers/multiphase/interFoam . 2) Rename the folder and the .C file. mv interFoam interTempFoam cd interTempFoam mv interFoam.C interTempFoam.C 3) Adapt the file Make/files to the changed naming and path. Code:
interTempFoam.C EXE = $(FOAM_USER_APPBIN)/interTempFoam In interTempFoam.C adapt the name of the following library being included. Code:
#include "myImmiscibleIncompressibleTwoPhaseMixture.H" Add the following lines in the file createFields.H. Code:
Info<< "Reading field U\n" << endl; volVectorField U ( IOobject ( "U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh ); //-------------modified_start--------------// Info<< "Reading field T\n" << endl; volScalarField T ( IOobject ( "T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh ); //-------------modified_end--------------// Code:
const dimensionedScalar& rho2 = mixture.rho2(); //-------------modified_start--------------// const dimensionedScalar& cp1 = mixture.cp1(); const dimensionedScalar& cp2 = mixture.cp2(); //not for Pr1 und Pr2 ??? //-------------modified_end--------------// Code:
// Mass flux surfaceScalarField rhoPhi ( IOobject ( "rhoPhi", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), fvc::interpolate(rho)*phi ); //-------------modified_start--------------// // Need to store rhoCp for ddt(rhoCp, T) volScalarField rhoCp ( IOobject ( "rhoCp", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT ), alpha1*rho1*cp1 + alpha2*rho2*cp2, alpha1.boundaryField().types() ); rhoCp.oldTime(); // heat flux surfaceScalarField rhoCpPhi ( IOobject ( "rhoCpPhi", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), fvc::interpolate(rhoCp)*phi ); //-------------modified_end--------------// In alphaEqn.H add: Code:
// Calculate the end-of-time-step mass flux rhoPhi = phiAlpha*(rho1 - rho2) + phi*rho2; //-------------modified_start--------------// // Calculate the end-of-time-step heat flux rhoCpPhi = phiAlpha*(rho1*cp1 - rho2*cp2) + phi*rho2*cp2; //-------------modified_end--------------// Code:
rho == alpha1*rho1 + alpha2*rho2; //-------------modified_start--------------// rhoCp == alpha1*rho1*cp1 + alpha2*rho2*cp2; //-------------modified_end--------------// Create a file named Teqn.H in the interTempFoam folder containing the Energy Equation: Code:
surfaceScalarField kappaf = incompressibleTwoPhaseMixture.kappaf(); fvScalarMatrix TEqn ( fvm::ddt(rhoCp,T) + fvm::div(rhoCpPhi,T) - fvm::laplacian(kappaf,T) ); TEqn.solve(); Add the command for solving the E-Eqn in interTempFoam.C. The specific location of the added line may be object of discussion. Here is my suggestion: Code:
// --- Pressure-velocity PIMPLE corrector loop while (pimple.loop()) { #include "alphaControls.H" #include "alphaEqnSubCycle.H" mixture.correct(); #include "UEqn.H" // --- Pressure corrector loop while (pimple.correct()) { #include "pEqn.H" } if (pimple.turbCorr()) { turbulence->correct(); } } //--------------modified_start----------------// #include "TEqn.H" //--------------modified_end------------------// Adapt the file Make/options to the names and locations of the beforehand modified libraries. The following code is probably not correct and thus object of my questions. Code:
EXE_INC = \ -I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \ -I$(LIB_SRC)/transportModels \ -I$(WM_PROJECT_USER_DIR)/src/transportModels/incompressible/lnInclude \ -I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \ -I$(LIB_SRC)/turbulenceModels/incompressible/turbulenceModel \ -I$(WM_PROJECT_USER_DIR)/src/transportModels/myImmiscibleIncompressibleTwoPhaseMixture/lnInclude \ -I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/fvOptions/lnInclude \ -I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/sampling/lnInclude EXE_LIBS = \ -L$(FOAM_USER_LIBBIN) \ -lmyImmiscibleIncompressibleTwoPhaseMixture \ -L$(FOAM_LIBBIN) \ -lincompressibleTurbulenceModel \ -lincompressibleRASModels \ -lincompressibleLESModels \ -lfiniteVolume \ -lfvOptions \ -lmeshTools \ -lsampling Compile the solver from within its directory interTempFoam. wclean wmake <-------------------------------------------------------------------------------------------------------------------------> When I do everything as described I get the following output log: Code:
~/OpenFOAM/puhlig-2.4.0/applications/solvers/multiphase/interTempFoam$ wmake Making dependency list for source file interTempFoam.C SOURCE=interTempFoam.C ; g++ -m64 -Dlinux64 -DWM_DP -Wall -Wextra -Wno-unused-parameter -Wold-style-cast -Wnon-virtual-dtor -O3 -DNoRepository -ftemplate-depth-100 -I/opt/openfoam240/src/transportModels/twoPhaseMixture/lnInclude -I/opt/openfoam240/src/transportModels -I/home/puhlig/OpenFOAM/puhlig-2.4.0/src/transportModels/incompressible/lnInclude -I/opt/openfoam240/src/transportModels/interfaceProperties/lnInclude -I/opt/openfoam240/src/turbulenceModels/incompressible/turbulenceModel -I/home/puhlig/OpenFOAM/puhlig-2.4.0/src/transportModels/myImmiscibleIncompressibleTwoPhaseMixture/lnInclude -I/opt/openfoam240/src/finiteVolume/lnInclude -I/opt/openfoam240/src/fvOptions/lnInclude -I/opt/openfoam240/src/meshTools/lnInclude -I/opt/openfoam240/src/sampling/lnInclude -IlnInclude -I. -I/opt/openfoam240/src/OpenFOAM/lnInclude -I/opt/openfoam240/src/OSspecific/POSIX/lnInclude -fPIC -c $SOURCE -o Make/linux64GccDPOpt/interTempFoam.o In file included from interTempFoam.C:108:0: TEqn.H: In function ‘int main(int, char**)’: TEqn.H:1:58: error: expected primary-expression before ‘.’ token surfaceScalarField kappaf = incompressibleTwoPhaseMixture.kappaf(); ^ In file included from interTempFoam.C:64:0: /opt/openfoam240/src/finiteVolume/lnInclude/readTimeControls.H:38:8: warning: unused variable ‘maxDeltaT’ [-Wunused-variable] scalar maxDeltaT = ^ make: *** [Make/linux64GccDPOpt/interTempFoam.o] Error 1 Do you have any idea how it comes to that error and how to solve it? I would very much appreciate comments from your side. Best Regards Perry Last edited by Pay; November 17, 2015 at 05:36. |
|
November 16, 2015, 20:57 |
interTempFoam in OF 2.3.1
|
#2 |
New Member
Hannah Dietterich
Join Date: Mar 2015
Posts: 8
Rep Power: 11 |
Hi,
I'm working on the exact same problem. I've used the wolfdynamics page to make nearly same modifications as you to incorporate temperature into OpenFOAM (2.3.1), and I am getting a similar error in the end. The only differences in my procedure are: Using "mixture" in your step 7: This is used to refer to the other quantities defined in myIncompressibleTwoPhaseMixture, so that's why I've used it. Anything else gives me the error you got. Code:
surfaceScalarField kappaf = mixture.kappaf(); fvScalarMatrix TEqn ( fvm::ddt(rhoCp, T) + fvm::div(rhoCpPhi, T) - fvm::laplacian(kappaf, T) ); TEqn.solve(); Code:
// --- Pressure-velocity PIMPLE corrector loop while (pimple.loop()) { #include "alphaControls.H" #include "alphaEqnSubCycle.H" mixture.correct(); #include "UEqn.H" #include "TEqn.H" // --- Pressure corrector loop while (pimple.correct()) { #include "pEqn.H" } if (pimple.turbCorr()) { turbulence->correct(); } } Code:
EXE_INC = \ -I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \ -I$(LIB_SRC)/transportModels \ -I$(WM_PROJECT_USER_DIR)/src/transportModels/incompressible/lnInclude \ -I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \ -I$(LIB_SRC)/turbulenceModels/incompressible/turbulenceModel \ -I$(WM_PROJECT_USER_DIR)/src/transportModels/myImmiscibleIncompressibleTwoPhaseMixture/lnInclude \ -I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/fvOptions/lnInclude \ -I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/sampling/lnInclude EXE_LIBS = \ -ltwoPhaseMixture \ -linterfaceProperties \ -ltwoPhaseProperties \ -L$(FOAM_USER_LIBBIN) \ -lmyIncompressibleTransportModels \ -lmyImmiscibleIncompressibleTwoPhaseMixture \ -lincompressibleTurbulenceModel \ -lincompressibleRASModels \ -lincompressibleLESModels \ -lfiniteVolume \ -lfvOptions \ -lmeshTools \ -lsampling Code:
In function `main': interThermalFoam.C:(.text.startup+0x409d): undefined reference to `Foam::incompressibleTwoPhaseMixture::kappaf() const' collect2: error: ld returned 1 exit status Does anyone have any suggestions? |
|
November 17, 2015, 03:00 |
|
#3 |
Senior Member
Anton Kidess
Join Date: May 2009
Location: Germany
Posts: 1,377
Rep Power: 30 |
Fantastic post Perry, you are a role model for this forum. Anyway, you did some strange things:
1. Your new library has the same name as the old one 2. You include both the new and old library in the Make files, although the old one is completely redundant.
__________________
*On twitter @akidTwit *Spend as much time formulating your questions as you expect people to spend on their answer. |
|
November 17, 2015, 11:04 |
|
#4 |
New Member
Perry Uhlig
Join Date: Oct 2015
Location: Germany
Posts: 13
Rep Power: 11 |
Hi Hannah,
thank you for your reply. I implemented the first line of your TEqn.H file: Code:
surfaceScalarField kappaf = mixture.kappaf(); I just adjusted the given damBreak tutorial case (with its very coarse mesh). It completed the calculation, but with temperature values beneath zero. In consequence I will try a refined mesh and I will probably have to look into fvSchemses and fvSolution in more detail. Greetings Perry |
|
November 17, 2015, 11:16 |
|
#5 | ||
New Member
Perry Uhlig
Join Date: Oct 2015
Location: Germany
Posts: 13
Rep Power: 11 |
Hello Anton,
thanks for your feedback and the compliment. Quote:
Code:
EXE_INC = \ -I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \ -I$(LIB_SRC)/transportModels \ -I$(WM_PROJECT_USER_DIR)/src/transportModels/incompressible/lnInclude \ - ... If NO: What do you mean? Quote:
Or didn't I? Greetings Perry |
|||
November 17, 2015, 12:50 |
|
#6 |
New Member
Hannah Dietterich
Join Date: Mar 2015
Posts: 8
Rep Power: 11 |
Hi Perry,
I'm glad that little change worked. I've modified mine to exactly match yours (with the mixture.kappaf) and it now compiles, too! Thanks for your very detailed description of all of your interFoam modifications. Now to see if it actually works... Hannah |
|
November 17, 2015, 13:28 |
|
#7 | |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
Quote:
Code:
fvm::ddt(T) + fvm::div(phi,T) + fvm::SuSp(-fvc::div(phi),T) - fvm::laplacian(alpha1*kappa1 / (rho1*Cp1), T) - fvm::laplacian((1 - alpha1)*kappa2 / (rho2*Cp2), T)
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
||
November 17, 2015, 15:20 |
|
#8 |
New Member
Hannah Dietterich
Join Date: Mar 2015
Posts: 8
Rep Power: 11 |
Unfortunately, I get a segmentation fault when I run my newly compiled interFoam with temperature for a modified damBreak case with temperature. How do I interpret the error message to identify what might be wrong?
Code:
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Create time Create mesh for time = 0 PIMPLE: Operating solver in PISO mode Reading field p_rgh Reading field U Reading field T Reading/calculating face flux field phi Reading transportProperties Selecting incompressible transport model Newtonian Selecting incompressible transport model Newtonian #0 Foam::error::printStack(Foam::Ostream&) at ??:? #1 Foam::sigSegv::sigHandler(int) at ??:? #2 at sigaction.c:? #3 at objectRegistry.C:? #4 Foam::objectRegistry::checkIn(Foam::regIOobject&) const at ??:? #5 Foam::regIOobject::checkIn() at ??:? #6 Foam::regIOobject::regIOobject(Foam::IOobject const&, bool) at ??:? #7 Foam::GeometricField<double, Foam::fvPatchField, Foam::volMesh>::GeometricField(Foam::IOobject const&, Foam::fvMesh const&, Foam::dimensioned<double> const&, Foam::word const&) at ??:? #8 Foam::interfaceProperties::interfaceProperties(Foam::GeometricField<double, Foam::fvPatchField, Foam::volMesh> const&, Foam::GeometricField<Foam::Vector<double>, Foam::fvPatchField, Foam::volMesh> const&, Foam::IOdictionary const&) at ??:? #9 Foam::immiscibleIncompressibleTwoPhaseMixture::immiscibleIncompressibleTwoPhaseMixture(Foam::GeometricField<Foam::Vector<double>, Foam::fvPatchField, Foam::volMesh> const&, Foam::GeometricField<double, Foam::fvsPatchField, Foam::surfaceMesh> const&) at ??:? #10 at ??:? #11 __libc_start_main at ??:? #12 at ??:? Segmentation fault Code:
div(rhoCpPhi,T) Gauss linearUpwind grad(T); div(rhoCpPhi) Gauss linearUpwind; Code:
T { solver BICCG; preconditioner DILU; tolerance 1e-07; relTol 0; } I'd appreciate help identifying where the error might be so that I might know what to do to fix it. mkraposhin, that's an interesting thought. Once this actually runs, then I'll be interested in getting the best solution possible. Hannah |
|
November 17, 2015, 16:29 |
|
#9 | |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
Quote:
immiscibleIncompressibleTwoPhaseMixture can you post your createFields.H file?
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
||
November 17, 2015, 16:31 |
|
#10 |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
Also, i think that it would be better to nest child class from mixture model, rather then creating new. Because at next lines OF creates turbulence model, that uses mixture model
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
|
November 17, 2015, 17:09 |
|
#11 | |
New Member
Hannah Dietterich
Join Date: Mar 2015
Posts: 8
Rep Power: 11 |
Quote:
Code:
Info<< "Reading field p_rgh\n" << endl; volScalarField p_rgh ( IOobject ( "p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh ); Info<< "Reading field U\n" << endl; volVectorField U ( IOobject ( "U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh ); // ADDED Info<< "Reading field T\n" << endl; volScalarField T ( IOobject ( "T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh ); // end ADDED #include "createPhi.H" Info<< "Reading transportProperties\n" << endl; immiscibleIncompressibleTwoPhaseMixture mixture(U, phi); volScalarField& alpha1(mixture.alpha1()); volScalarField& alpha2(mixture.alpha2()); const dimensionedScalar& rho1 = mixture.rho1(); const dimensionedScalar& rho2 = mixture.rho2(); // ADDED const dimensionedScalar& cp1 = mixture.cp1(); const dimensionedScalar& cp2 = mixture.cp2(); //const dimensionedScalar& Pr1 = mixture.Pr1(); //const dimensionedScalar& Pr2 = mixture.Pr2(); // end ADDED // Need to store rho for ddt(rho, U) volScalarField rho ( IOobject ( "rho", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT ), alpha1*rho1 + alpha2*rho2, alpha1.boundaryField().types() ); rho.oldTime(); // ADDED // Need to store rhoCp volScalarField rhoCp ( IOobject ( "rhoCp", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT ), alpha1*rho1*cp1 + alpha2*rho2*cp2, alpha1.boundaryField().types() ); rhoCp.oldTime(); // end ADDED // Mass flux surfaceScalarField rhoPhi ( IOobject ( "rhoPhi", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), fvc::interpolate(rho)*phi ); // ADDED surfaceScalarField rhoCpPhi ( IOobject ( "rhoCpPhi", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), fvc::interpolate(rhoCp)*phi ); // end ADDED // Construct incompressible turbulence model autoPtr<incompressible::turbulenceModel> turbulence ( incompressible::turbulenceModel::New(U, phi, mixture) ); #include "readGravitationalAcceleration.H" /* dimensionedVector g0(g); // Read the data file and initialise the interpolation table interpolationTable<vector> timeSeriesAcceleration ( runTime.path()/runTime.caseConstant()/"acceleration.dat" ); */ Info<< "Calculating field g.h\n" << endl; volScalarField gh("gh", g & mesh.C()); surfaceScalarField ghf("ghf", g & mesh.Cf()); volScalarField p ( IOobject ( "p", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE ), p_rgh + rho*gh ); label pRefCell = 0; scalar pRefValue = 0.0; setRefCell ( p, p_rgh, mesh.solutionDict().subDict("PIMPLE"), pRefCell, pRefValue ); if (p_rgh.needReference()) { p += dimensionedScalar ( "p", p.dimensions(), pRefValue - getRefCellValue(p, pRefCell) ); p_rgh = p - rho*gh; } fv::IOoptionList fvOptions(mesh); // Creating field for thermophysical properties // MULES Correction tmp<surfaceScalarField> tphiAlphaCorr0; Code:
#include "myIncompressibleTwoPhaseMixture.H" #include "addToRunTimeSelectionTable.H" #include "surfaceFields.H" #include "fvc.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(incompressibleTwoPhaseMixture, 0); } |
||
November 18, 2015, 01:22 |
|
#12 | |
Senior Member
Anton Kidess
Join Date: May 2009
Location: Germany
Posts: 1,377
Rep Power: 30 |
Sorry for the confusion, I just realized I didn't read it properly and indeed you don't have the original immiscibleIncompressibleTwoPhaseMixture in Make/*.
Quote:
__________________
*On twitter @akidTwit *Spend as much time formulating your questions as you expect people to spend on their answer. |
||
November 18, 2015, 02:22 |
|
#13 | |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
Quote:
So, you need to change immiscibleIncompressibleTwoPhaseMixture to myImmiscibleIncompressibleTwoPhaseMixture But when you do that, i think, that solver will fail to compile, because turbulence model class depends on mixture class. Later i can try to explain how to derive your own class from class immiscibleIncompressibleTwoPhaseMixture
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
||
November 18, 2015, 02:44 |
|
#14 |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
O.k., i'm sorry,
i didn't read untill the end of your post, but nevertheless, turbulence model needs original data type, but not modified, that's why it fails to run and of course, you must avoid using same data type names for different classes
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
|
November 18, 2015, 06:46 |
|
#15 |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
Hi, i spend ~ 20 minutes and made child class myImmiscibleIncompressibleTwoPhaseMixture, which can do all work without substituting user type instead of original OF
Step 1. Create your new class myImmiscibleIncompressibleTwoPhaseMixture as a sub class of immiscibleIncompressibleTwoPhaseMixture, see attachment files .H and .C. Please note that for compatibility with OF you need to change all names, for example: Code:
#ifndef myImmiscibleIncompressibleTwoPhaseMixture_H #define myImmiscibleIncompressibleTwoPhaseMixture_H Step 2. Adjust files Make/files and Make/options to include all neccesary files from OF Make/files Code:
myImmiscibleIncompressibleTwoPhaseMixture.C LIB = $(FOAM_USER_LIBBIN)/libmyImmiscibleIncompressibleTwoPhaseMixture Code:
EXE_INC = \ -I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/OpenFOAM/lnInclude \ -I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/transportModels \ -I$(LIB_SRC)/transportModels/immiscibleIncompressibleTwoPhaseMixture/lnInclude \ -I$(LIB_SRC)/transportModels/incompressible/lnInclude \ -I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \ -I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \ -I$(LIB_SRC)/finiteVolume/lnInclude LIB_LIBS = \ -ltwoPhaseMixture \ -lincompressibleTransportModels \ -linterfaceProperties \ -ltwoPhaseProperties \ -limmiscibleIncompressibleTwoPhaseMixture \ -lfiniteVolume Step 4. Initialize new model in new solver: - add -I../myImmiscibleIncompressibleTwoPhaseMixture \ line to Make/options EXE_INC section of your solver - add consequent lines to -L section of Make/options file of your solver - add #include "myImmiscibleIncompressibleTwoPhaseMixture.H" to solver instead of OF original solver - add initialization of your myImmiscibleIncompressibleTwoPhaseMixture class instead of OF original class in createFields.H - use method mixture.lambdaEff() to get mixture heat conduction coefficient - solve equation for temperature, formulated in volumetric fluxes: Code:
solve ( fvm::ddt(T) + fvm::div(phi,T) + fvm::SuSp(-fvc::div(phi),T) - fvm::laplacian(lambdam, T) )
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin Last edited by mkraposhin; November 19, 2015 at 13:00. |
|
November 18, 2015, 13:07 |
|
#16 |
New Member
Hannah Dietterich
Join Date: Mar 2015
Posts: 8
Rep Power: 11 |
Thanks so much for your efforts. I'll go through and change all the internal names as well. Since all the new properties are defined in myIncompressibleTwoPhaseMixture (the immiscible one just calls them), that will need to be adjusted, too.
Edit: I see you've suggested moving all of the definitions into myImmiscibleIncompressibleTwoPhaseMixture. I'll incorporate that. This doesn't help explain, though, why this isn't needed for Perry and everyone else using these (http://www.wolfdynamics.com/images/c...erTempFoam.pdf, http://www.cfd-online.com/Forums/ope...interfoam.html, etc) instructions. Hannah |
|
November 18, 2015, 18:15 |
|
#17 |
New Member
Hannah Dietterich
Join Date: Mar 2015
Posts: 8
Rep Power: 11 |
Well those changes just produced a lot of compiling errors, but continuing to mess around with my original Make/options and #include's produced a version that compiled and ran a test case successfully. I'll look into making mkraposhin's suggested modifications to TEqn.
In the meantime, I'm attaching all of the files in case they might help someone else (using OF 2.3.1). Hannah |
|
November 19, 2015, 05:27 |
|
#18 |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
I will try to compile your code and report what is wrong later.
Regarding previous version of library and solver - i made it for OF 2.4.0, not for 2.3.1
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
|
November 19, 2015, 06:28 |
|
#19 | ||||
New Member
Perry Uhlig
Join Date: Oct 2015
Location: Germany
Posts: 13
Rep Power: 11 |
Hello Matvey,
first of all, thanks for taking part in the thread. Quote:
in alphaEqn.H: Code:
// Calculate the end-of-time-step mass flux rhoPhi = phiAlpha*(rho1 - rho2) + phi*rho2; } //-------------modified_start--------------// // Calculate the end-of-time-step heat flux rhoCpPhi = phiAlpha*(rho1*cp1 - rho2*cp2) + phi*rho2*cp2; //-------------modified_end--------------// Code:
rho == alpha1*rho1 + alpha2*rho2; //-------------modified_start--------------// rhoCp == alpha1*rho1*cp1 + alpha2*rho2*cp2; //-------------modified_end--------------// I see that you divided the equation by (rho*cp) viewing the terms separately: Quote:
Quote:
Quote:
What does this new term do? I do not know the operation in there (SuSp(-fvc ...). Best Regards Perry |
|||||
November 19, 2015, 13:35 |
|
#20 |
Senior Member
Matvey Kraposhin
Join Date: Mar 2009
Location: Moscow, Russian Federation
Posts: 355
Rep Power: 21 |
@ Hannah:
Hi, i looked at your code and i understand that you mixed from all versions So, i used my approach to create a new class for mixture and mass flux approach for energy equation. It works for OF 2.3.1. You can find source code in the attachment. To make solver, go into folder interThermalFoam/myImmiscibleIncompressibleTwoPhaseMixture and run wmake libso then go one folder up and run wmake
__________________
MDPI Fluids (Q2) special issue for OSS software: https://www.mdpi.com/journal/fluids/..._modelling_OSS GitHub: https://github.com/unicfdlab Linkedin: https://linkedin.com/in/matvey-kraposhin-413869163 RG: https://www.researchgate.net/profile/Matvey_Kraposhin |
|
Tags |
2.4.0, compile problems, interfoam, teqn |
Thread Tools | Search this Thread |
Display Modes | |
|
|
Similar Threads | ||||
Thread | Thread Starter | Forum | Replies | Last Post |
Source Term due to evaporation in energy transport equation | styleworker | OpenFOAM Programming & Development | 3 | September 7, 2022 04:09 |
Adding Energy Equation for porousSimpleFoam | yesaswi92 | OpenFOAM Programming & Development | 20 | April 13, 2018 09:17 |
Definition of energy equation using UDF | aestas | Fluent UDF and Scheme Programming | 15 | April 28, 2015 04:33 |
help adding the energy equation to porousinterfoam using the enthalpy formulation | nadine | OpenFOAM Programming & Development | 18 | June 17, 2014 09:39 |
Why FVM for high-Re flows? | Zhong Lei | Main CFD Forum | 23 | May 14, 1999 14:22 |