View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0002323 | OpenFOAM | Bug | public | 2016-11-08 15:00 | 2018-07-10 11:22 |
Reporter | nima.samkhaniani | Assigned To | henry | ||
Priority | normal | Severity | crash | Reproducibility | have not tried |
Status | resolved | Resolution | fixed | ||
Platform | GNU/Linux | OS | Ubuntu | OS Version | 12.04 |
Fixed in Version | dev | ||||
Summary | 0002323: thermalBaffel1D | ||||
Description | I run a dynamicMesh problem similar movingCone using compressibelInterDyMFoam. I noticed that it failed after several iterations when i uses thermalBaffel1D, I dont know why but the size of scalarField kappas and baffleThickness() in thermalBaffel1D becomes unequal and it is the source of error. i solved this problem via replacing The following code in thermalBaffle1DFvPatchScalarField.C KDeltaSolid = kappas/baffleThickness(); with this scalarField KDeltaSolid(patch().size(), 0.0); scalarField alpha(patch().size(), 0.0); scalarField baffleThickness1=baffleThickness(); forAll(kappas, i) { kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0); KDeltaSolid[i] = kappas[i]/baffleThickness1[i]; //add } | ||||
Tags | No tags attached. | ||||
|
1) Can you attach a testcase? 2) Have you run your case through valgrind? |
|
1)The file is attached. It crashes before time 0.038 2)nope |
|
The problem is the underlying 'mapped' geometric matching which does not support topology changes out of the box so the mapping would still try to get data for 20 faces even though the patch was now 19 faces. Attached fix for (dev) src/TurbulenceModels/compressible/turbulentFluidThermoModels/derivedFvPatchFields/thermalBaffle1D/thermalBaffle1DFvPatchScalarField.C will get around this. We will have to think about where to fix this properly. thermalBaffle1DFvPatchScalarField.C (12,230 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2016 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 "volFields.H" #include "surfaceFields.H" #include "mappedPatchBase.H" #include "turbulentFluidThermoModel.H" #include "mapDistribute.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { namespace compressible { // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template<class solidType> thermalBaffle1DFvPatchScalarField<solidType>:: thermalBaffle1DFvPatchScalarField ( const fvPatch& p, const DimensionedField<scalar, volMesh>& iF ) : mappedPatchBase(p.patch()), mixedFvPatchScalarField(p, iF), TName_("T"), baffleActivated_(true), thickness_(p.size()), Qs_(p.size()), solidDict_(), solidPtr_(nullptr), QrPrevious_(p.size()), QrRelaxation_(1), QrName_("undefined-Qr") {} template<class solidType> thermalBaffle1DFvPatchScalarField<solidType>:: thermalBaffle1DFvPatchScalarField ( const thermalBaffle1DFvPatchScalarField& ptf, const fvPatch& p, const DimensionedField<scalar, volMesh>& iF, const fvPatchFieldMapper& mapper ) : mappedPatchBase(p.patch(), ptf), mixedFvPatchScalarField(ptf, p, iF, mapper), TName_(ptf.TName_), baffleActivated_(ptf.baffleActivated_), thickness_(ptf.thickness_, mapper), Qs_(ptf.Qs_, mapper), solidDict_(ptf.solidDict_), solidPtr_(ptf.solidPtr_), QrPrevious_(ptf.QrPrevious_, mapper), QrRelaxation_(ptf.QrRelaxation_), QrName_(ptf.QrName_) {} template<class solidType> thermalBaffle1DFvPatchScalarField<solidType>:: thermalBaffle1DFvPatchScalarField ( const fvPatch& p, const DimensionedField<scalar, volMesh>& iF, const dictionary& dict ) : mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict), mixedFvPatchScalarField(p, iF), TName_("T"), baffleActivated_(dict.lookupOrDefault<bool>("baffleActivated", true)), thickness_(), Qs_(p.size(), 0), solidDict_(dict), solidPtr_(), QrPrevious_(p.size(), 0.0), QrRelaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)), QrName_(dict.lookupOrDefault<word>("Qr", "none")) { fvPatchScalarField::operator=(scalarField("value", dict, p.size())); if (dict.found("thickness")) { thickness_ = scalarField("thickness", dict, p.size()); } if (dict.found("Qs")) { Qs_ = scalarField("Qs", dict, p.size()); } if (dict.found("QrPrevious")) { QrPrevious_ = scalarField("QrPrevious", dict, p.size()); } if (dict.found("refValue") && baffleActivated_) { // Full restart refValue() = scalarField("refValue", dict, p.size()); refGrad() = scalarField("refGradient", dict, p.size()); valueFraction() = scalarField("valueFraction", dict, p.size()); } else { // Start from user entered data. Assume zeroGradient. refValue() = *this; refGrad() = 0.0; valueFraction() = 0.0; } } template<class solidType> thermalBaffle1DFvPatchScalarField<solidType>:: thermalBaffle1DFvPatchScalarField ( const thermalBaffle1DFvPatchScalarField& ptf ) : mappedPatchBase(ptf.patch().patch(), ptf), mixedFvPatchScalarField(ptf), TName_(ptf.TName_), baffleActivated_(ptf.baffleActivated_), thickness_(ptf.thickness_), Qs_(ptf.Qs_), solidDict_(ptf.solidDict_), solidPtr_(ptf.solidPtr_), QrPrevious_(ptf.QrPrevious_), QrRelaxation_(ptf.QrRelaxation_), QrName_(ptf.QrName_) {} template<class solidType> thermalBaffle1DFvPatchScalarField<solidType>:: thermalBaffle1DFvPatchScalarField ( const thermalBaffle1DFvPatchScalarField& ptf, const DimensionedField<scalar, volMesh>& iF ) : mappedPatchBase(ptf.patch().patch(), ptf), mixedFvPatchScalarField(ptf, iF), TName_(ptf.TName_), baffleActivated_(ptf.baffleActivated_), thickness_(ptf.thickness_), Qs_(ptf.Qs_), solidDict_(ptf.solidDict_), solidPtr_(ptf.solidPtr_), QrPrevious_(ptf.QrPrevious_), QrRelaxation_(ptf.QrRelaxation_), QrName_(ptf.QrName_) {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template<class solidType> bool thermalBaffle1DFvPatchScalarField<solidType>::owner() const { const label patchi = patch().index(); const label nbrPatchi = samplePolyPatch().index(); return (patchi < nbrPatchi); } template<class solidType> const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid() const { if (this->owner()) { if (solidPtr_.empty()) { solidPtr_.reset(new solidType(solidDict_)); } return solidPtr_(); } else { const fvPatch& nbrPatch = patch().boundaryMesh()[samplePolyPatch().index()]; const thermalBaffle1DFvPatchScalarField& nbrField = refCast<const thermalBaffle1DFvPatchScalarField> ( nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_) ); return nbrField.solid(); } } template<class solidType> tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>:: baffleThickness() const { if (this->owner()) { if (thickness_.size() != patch().size()) { FatalIOErrorInFunction ( solidDict_ )<< " Field thickness has not been specified " << " for patch " << this->patch().name() << exit(FatalIOError); } return thickness_; } else { const mapDistribute& mapDist = this->mappedPatchBase::map(); const fvPatch& nbrPatch = patch().boundaryMesh()[samplePolyPatch().index()]; const thermalBaffle1DFvPatchScalarField& nbrField = refCast<const thermalBaffle1DFvPatchScalarField> ( nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_) ); tmp<scalarField> tthickness ( new scalarField(nbrField.baffleThickness()) ); scalarField& thickness = tthickness.ref(); mapDist.distribute(thickness); return tthickness; } } template<class solidType> tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::Qs() const { if (this->owner()) { return Qs_; } else { const mapDistribute& mapDist = this->mappedPatchBase::map(); const fvPatch& nbrPatch = patch().boundaryMesh()[samplePolyPatch().index()]; const thermalBaffle1DFvPatchScalarField& nbrField = refCast<const thermalBaffle1DFvPatchScalarField> ( nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_) ); tmp<scalarField> tQs(new scalarField(nbrField.Qs())); scalarField& Qs = tQs.ref(); mapDist.distribute(Qs); return tQs; } } template<class solidType> void thermalBaffle1DFvPatchScalarField<solidType>::autoMap ( const fvPatchFieldMapper& m ) { // Clear out addressing in case of topo change this->mappedPatchBase::clearOut(); mixedFvPatchScalarField::autoMap(m); if (this->owner()) { thickness_.autoMap(m); Qs_.autoMap(m); } } template<class solidType> void thermalBaffle1DFvPatchScalarField<solidType>::rmap ( const fvPatchScalarField& ptf, const labelList& addr ) { mixedFvPatchScalarField::rmap(ptf, addr); const thermalBaffle1DFvPatchScalarField& tiptf = refCast<const thermalBaffle1DFvPatchScalarField>(ptf); if (this->owner()) { thickness_.rmap(tiptf.thickness_, addr); Qs_.rmap(tiptf.Qs_, addr); } } template<class solidType> void thermalBaffle1DFvPatchScalarField<solidType>::updateCoeffs() { if (updated()) { return; } // Since we're inside initEvaluate/evaluate there might be processor // comms underway. Change the tag we use. int oldTag = UPstream::msgType(); UPstream::msgType() = oldTag+1; const mapDistribute& mapDist = this->mappedPatchBase::map(); const label patchi = patch().index(); const label nbrPatchi = samplePolyPatch().index(); if (baffleActivated_) { const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi]; const compressible::turbulenceModel& turbModel = db().template lookupObject<compressible::turbulenceModel> ( turbulenceModel::propertiesName ); // local properties const scalarField kappaw(turbModel.kappaEff(patchi)); const fvPatchScalarField& Tp = patch().template lookupPatchField<volScalarField, scalar>(TName_); scalarField Qr(Tp.size(), 0.0); if (QrName_ != "none") { Qr = patch().template lookupPatchField<volScalarField, scalar> (QrName_); Qr = QrRelaxation_*Qr + (1.0 - QrRelaxation_)*QrPrevious_; QrPrevious_ = Qr; } tmp<scalarField> Ti = patchInternalField(); scalarField myKDelta(patch().deltaCoeffs()*kappaw); // nrb properties scalarField nbrTp = turbModel.transport().T().boundaryField()[nbrPatchi]; mapDist.distribute(nbrTp); // solid properties scalarField kappas(patch().size(), 0.0); forAll(kappas, i) { kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0); } scalarField KDeltaSolid(kappas/baffleThickness()); scalarField alpha(KDeltaSolid - Qr/Tp); valueFraction() = alpha/(alpha + myKDelta); refValue() = (KDeltaSolid*nbrTp + Qs()/2.0)/alpha; if (debug) { scalar Q = gAverage(kappaw*snGrad()); Info<< patch().boundaryMesh().mesh().name() << ':' << patch().name() << ':' << this->internalField().name() << " <- " << nbrPatch.name() << ':' << this->internalField().name() << " :" << " heat[W]:" << Q << " walltemperature " << " min:" << gMin(*this) << " max:" << gMax(*this) << " avg:" << gAverage(*this) << endl; } } // Restore tag UPstream::msgType() = oldTag; mixedFvPatchScalarField::updateCoeffs(); } template<class solidType> void thermalBaffle1DFvPatchScalarField<solidType>::write(Ostream& os) const { mixedFvPatchScalarField::write(os); mappedPatchBase::write(os); if (this->owner()) { baffleThickness()().writeEntry("thickness", os); Qs()().writeEntry("Qs", os); solid().write(os); } QrPrevious_.writeEntry("QrPrevious", os); os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl; os.writeKeyword("relaxation")<< QrRelaxation_ << token::END_STATEMENT << nl; } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace compressible } // End namespace Foam // ************************************************************************* // |
|
Resolved in OpenFOAM-dev by commit 4a45e4f519feef615630ca50dc13f8f0517c17f9 |
Date Modified | Username | Field | Change |
---|---|---|---|
2016-11-08 15:00 | nima.samkhaniani | New Issue | |
2016-11-08 15:00 | nima.samkhaniani | Tag Attached: thermal baffel + dynamicMesh | |
2016-11-11 15:23 | MattijsJ | Note Added: 0007151 | |
2016-11-21 18:22 | nima.samkhaniani | File Added: movingPistonBaffel-compressibleInterDyMFoam.zip | |
2016-11-21 18:22 | nima.samkhaniani | Note Added: 0007270 | |
2016-11-22 15:02 | MattijsJ | File Added: thermalBaffle1DFvPatchScalarField.C | |
2016-11-22 15:02 | MattijsJ | Note Added: 0007275 | |
2016-12-05 08:49 | henry | Assigned To | => henry |
2016-12-05 08:49 | henry | Status | new => resolved |
2016-12-05 08:49 | henry | Resolution | open => fixed |
2016-12-05 08:49 | henry | Fixed in Version | => dev |
2016-12-05 08:49 | henry | Note Added: 0007404 | |
2018-07-10 11:22 | administrator | Tag Detached: thermal baffel + dynamicMesh |