View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0001902 | OpenFOAM | Bug | public | 2015-11-06 08:21 | 2016-02-01 13:05 |
Reporter | projectionist | Assigned To | henry | ||
Priority | normal | Severity | minor | Reproducibility | always |
Status | resolved | Resolution | fixed | ||
Platform | GNU/Linux | OS | Ubuntu | OS Version | 15.04 |
Product Version | dev | ||||
Summary | 0001902: [reactingTwoPhaseEulerFoam] pureIsothermalPhaseModel not usable | ||||
Description | When selecting the pureIsothermalPhaseModel, I get the following error message: --> FOAM FATAL ERROR: Not implemented From function virtual const Foam::tmp<Foam::GeometricField<double, Foam::fvPatchField, Foam::volMesh> >& Foam::phaseModel::divU() const in file phaseModel/phaseModel/phaseModel.C at line 170. FOAM aborting #0 Foam::error::printStack(Foam::Ostream&) at ??:? #1 Foam::error::abort() at ??:? #2 Foam::phaseModel::divU() const at ??:? #3 Foam::twoPhaseSystem::solve() at ??:? #4 ? at ??:? #5 __libc_start_main in "/lib/x86_64-linux-gnu/libc.so.6" #6 ? at ??:? Abgebrochen (Speicherabzug geschrieben) | ||||
Steps To Reproduce | When twoPhaseSystem.solve() is called, the divU() method of the phaseModel gets called. Since, pureIsothermalPhaseModel does not overwrite the (I guess placeholder-) implementation of divU() of the base class phaseModel, it cannot be used. const Foam::tmp<Foam::volScalarField>& Foam::phaseModel::divU() const { NotImplemented; static tmp<Foam::volScalarField> divU_(NULL); return divU_; } | ||||
Additional Information | The responsible code of solve() from twoPhaseSystem.C is here: if (phase1_.divU().valid() && phase2_.divU().valid()) { tdgdt = ( alpha2.dimensionedInternalField() *phase1_.divU()().dimensionedInternalField() - alpha1.dimensionedInternalField() *phase2_.divU()().dimensionedInternalField() ); } else if (phase1_.divU().valid()) { tdgdt = ( alpha2.dimensionedInternalField() *phase1_.divU()().dimensionedInternalField() ); } else if (phase2_.divU().valid()) { tdgdt = ( - alpha1.dimensionedInternalField() *phase2_.divU()().dimensionedInternalField() ); } The ad-hoc solution might be to get pureIsothermalPhaseModel.divU() to return an invalid volScalarField. | ||||
Tags | No tags attached. | ||||
|
Do you need the IsothermalPhaseModel to support compressibility? This may be useful is the case of a bubble column in which the bubbles expand as they rise but assumed to be in thermal equilibrium with the isothermal liquid. In order to specifically handle the case of incompressible, isothermal phases we could implement IsochoricPhaseModel if that is also needed. |
|
IsothermalPhaseModel.H (2,980 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2015 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 <http://www.gnu.org/licenses/>. Class Foam::IsothermalPhaseModel Description Class which represents a phase for which the temperature (strictly energy) remains constant. Returns an empty energy equation and does nothing when correctThermo is called. SourceFiles IsothermalPhaseModel.C \*---------------------------------------------------------------------------*/ #ifndef IsothermalPhaseModel_H #define IsothermalPhaseModel_H #include "phaseModel.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { /*---------------------------------------------------------------------------*\ Class phaseModel Declaration \*---------------------------------------------------------------------------*/ template<class BasePhaseModel> class IsothermalPhaseModel : public BasePhaseModel { // Private data //- Dilatation rate tmp<volScalarField> divU_; public: // Constructors IsothermalPhaseModel ( const phaseSystem& fluid, const word& phaseName, const label index ); //- Destructor virtual ~IsothermalPhaseModel(); // Member Functions //- Correct the thermodynamics virtual void correctThermo(); //- Return the enthalpy equation virtual tmp<fvScalarMatrix> heEqn(); //- Return an invalid tmp-pointer virtual const tmp<volScalarField>& divU() const; }; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace Foam // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #ifdef NoRepository # include "IsothermalPhaseModel.C" #endif // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #endif // ************************************************************************* // |
|
IsothermalPhaseModel.C (2,294 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2015 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 <http://www.gnu.org/licenses/>. \*---------------------------------------------------------------------------*/ #include "IsothermalPhaseModel.H" #include "phaseSystem.H" // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template<class BasePhaseModel> Foam::IsothermalPhaseModel<BasePhaseModel>::IsothermalPhaseModel ( const phaseSystem& fluid, const word& phaseName, const label index ) : BasePhaseModel(fluid, phaseName, index), divU_(NULL) {} // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // template<class BasePhaseModel> Foam::IsothermalPhaseModel<BasePhaseModel>::~IsothermalPhaseModel() {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template<class BasePhaseModel> void Foam::IsothermalPhaseModel<BasePhaseModel>::correctThermo() {} template<class BasePhaseModel> Foam::tmp<Foam::fvScalarMatrix> Foam::IsothermalPhaseModel<BasePhaseModel>::heEqn() { return tmp<fvScalarMatrix>(); } template<class BasePhaseModel> const Foam::tmp<Foam::volScalarField>& Foam::IsothermalPhaseModel<BasePhaseModel>::divU() const { return divU_; } // ************************************************************************* // |
|
The attached files do the trick for me: In the header file I added a new private data member (tmp<volScalarField> divU_), which is initialized by the constructor (divU_(NULL)). Furthermore, the method divU() was added: virtual const tmp<volScalarField>& divU() const; In the implementation file, this additional method just passes a reference to divU_: template<class BasePhaseModel> const Foam::tmp<Foam::volScalarField>& Foam::IsothermalPhaseModel<BasePhaseModel>::divU() const { return divU_; } |
|
I intend to use this phaseModel to suppress the solution of the energy equation, which causes me headaches on covergence or even crashes my simulations via the "Maximum number of iterations exceeded" error. Up to now, my ad-hoc solution to iso-thermal problems, was to clone the solver and comment the lines which include the EEqns.H file. A very inelegant solution, but it worked. I am under the impression, that the IsothermalPhaseModel is the solution to my problem. Since, the solvers are compressible by design, it is the phaseModel that introduces incompressibility again by providing an empty energy equation. However, I could be wrong in my assessment. |
|
Sorry for again not answering your question. The IsochoricPhaseModel is indeed the solution that would fit my problem. |
|
The problem with the IsochoricPhaseModel concept is that it needn't be isothermal and IsothermalPhaseModel needn't be isochoric. In your simulations do you need the option to include the variation bubble-size with pressure? If you could provide a bit more detail about your cases it would make it easier to work out the best approach. |
|
Right now, I am simulating a small bubble column at room temperature. Thus, heat transfer and bubble size change are neglected. For my larger cases bubble size variation due to pressure, might be useful. So I guess, heat transfer can be more readily neglected than buble size change. |
|
I've noticed that it is often useful to suppress the solution of the energy equations when the simulation is initialized. When the the initial pressure shocks have dissipated the energy solution can be turned back on. In my modified solver I have an extra corrector loop over the energy equations, so the energy solution can simply be temporarily disabled by setting nEnergyCorr 0. |
|
I have moved the dilatation into the MovingPhaseModel to support phase volume fraction changes due to pressure: commit fdc9432c4a3be1b0b833ebcf423e1b541801f09c This should also resolve the problem you are having with the IsothermalPhaseModel. |
|
resolved by commit fdc9432c4a3be1b0b833ebcf423e1b541801f09c |
Date Modified | Username | Field | Change |
---|---|---|---|
2015-11-06 08:21 | projectionist | New Issue | |
2015-11-06 08:36 | henry | Note Added: 0005567 | |
2015-11-06 08:44 | projectionist | File Added: IsothermalPhaseModel.H | |
2015-11-06 08:44 | projectionist | File Added: IsothermalPhaseModel.C | |
2015-11-06 08:47 | projectionist | Note Added: 0005568 | |
2015-11-06 08:54 | projectionist | Note Added: 0005569 | |
2015-11-06 08:58 | projectionist | Note Added: 0005570 | |
2015-11-06 09:39 | henry | Note Added: 0005571 | |
2015-11-06 10:05 | projectionist | Note Added: 0005572 | |
2015-11-06 11:55 | Juho | Note Added: 0005573 | |
2015-11-06 15:39 | henry | Note Added: 0005576 | |
2015-11-10 09:12 | henry | Note Added: 0005596 | |
2015-11-10 09:12 | henry | Status | new => resolved |
2015-11-10 09:12 | henry | Resolution | open => fixed |
2015-11-10 09:12 | henry | Assigned To | => henry |