Loading...
Searching...
No Matches
CoulombForce.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2023-2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "CoulombForce.H"
29#include "demandDrivenData.H"
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class CloudType>
34Foam::volVectorField& Foam::CoulombForce<CloudType>::getOrReadField
35(
36 const word& fieldName
37) const
38{
39 auto* ptr = this->mesh().template getObjectPtr<volVectorField>(fieldName);
40
41 if (!ptr)
42 {
43 ptr = new volVectorField
44 (
45 IOobject
46 (
47 fieldName,
48 this->mesh().time().timeName(),
49 this->mesh().thisDb(),
50 IOobject::MUST_READ,
51 IOobject::AUTO_WRITE,
52 IOobject::REGISTER
53 ),
54 this->mesh()
55 );
56 regIOobject::store(ptr);
57 }
58
59 return *ptr;
60}
61
62
63// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64
65template<class CloudType>
67(
68 CloudType& owner,
69 const fvMesh& mesh,
70 const dictionary& dict
71)
72:
73 ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
74 qPtr_
75 (
76 Function1<scalar>::New("q", this->coeffs(), &mesh)
77 ),
78 Ename_(this->coeffs().template getOrDefault<word>("E", "E")),
79 EInterpPtr_(nullptr)
80{}
81
82
83template<class CloudType>
85(
86 const CoulombForce& pf
87)
88:
90 qPtr_(pf.qPtr_.clone()),
91 Ename_(pf.Ename_),
92 EInterpPtr_(nullptr)
93{}
94
95
96// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97
98template<class CloudType>
100{
101 if (store)
102 {
103 const volVectorField& E = getOrReadField(Ename_);
104
105 EInterpPtr_.reset
106 (
108 (
109 this->owner().solution().interpolationSchemes(),
110 E
111 ).ptr()
112 );
113 }
114 else
116 EInterpPtr_.reset(nullptr);
117 }
118}
119
120
121template<class CloudType>
123(
124 const typename CloudType::parcelType& p,
125 const typename CloudType::parcelType::trackingData& td,
126 const scalar dt,
127 const scalar mass,
128 const scalar Re,
129 const scalar muc
130) const
131{
132 const interpolation<vector>& EInterp = *EInterpPtr_;
133
134 const scalar q = qPtr_->value(p.d());
135
136 // (YSSD:Eq. 6 - left term)
137 return forceSuSp
138 (
139 q*EInterp.interpolate(p.coordinates(), p.currentTetIndices()),
140 Zero
141 );
142}
143
144
145// ************************************************************************* //
Particle electric force model involving the Coulomb force calculation.
virtual autoPtr< ParticleForce< CloudType > > clone() const
Construct and return a clone.
virtual forceSuSp calcNonCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
virtual void cacheFields(const bool store)
Cache fields.
CoulombForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
Abstract base class for particle forces.
const fvMesh & mesh() const noexcept
Return the mesh database.
static autoPtr< ParticleForce< CloudType > > New(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType)
Selector.
const CloudType & owner() const noexcept
Return const access to the cloud owner.
ParticleForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType, const bool readCoeffs)
Construct from mesh.
const dictionary & coeffs() const noexcept
Return the force coefficients dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Helper container for force Su and Sp terms.
Definition forceSuSp.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for volume field interpolation.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
word timeName
Definition getTimeIndex.H:3
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dictionary dict