Loading...
Searching...
No Matches
ParamagneticForce.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "demandDrivenData.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class CloudType>
37(
39 const fvMesh& mesh,
40 const dictionary& dict
41)
42:
44 HdotGradHName_
45 (
46 this->coeffs().template getOrDefault<word>("HdotGradH", "HdotGradH")
47 ),
48 HdotGradHInterpPtr_(nullptr),
49 magneticSusceptibility_
50 (
51 this->coeffs().getScalar("magneticSusceptibility")
52 )
53{}
54
55
56template<class CloudType>
58(
59 const ParamagneticForce& pf
60)
61:
63 HdotGradHName_(pf.HdotGradHName_),
64 HdotGradHInterpPtr_(pf.HdotGradHInterpPtr_),
65 magneticSusceptibility_(pf.magneticSusceptibility_)
66{}
67
68
69// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
70
71template<class CloudType>
73{}
74
75
76// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77
78template<class CloudType>
80{
81 if (store)
82 {
83 const volVectorField& HdotGradH =
84 this->mesh().template lookupObject<volVectorField>(HdotGradHName_);
85
86 HdotGradHInterpPtr_ = interpolation<vector>::New
87 (
88 this->owner().solution().interpolationSchemes(),
89 HdotGradH
90 ).ptr();
91 }
92 else
93 {
94 deleteDemandDrivenData(HdotGradHInterpPtr_);
95 }
96}
97
98
99template<class CloudType>
101(
102 const typename CloudType::parcelType& p,
103 const typename CloudType::parcelType::trackingData& td,
104 const scalar dt,
105 const scalar mass,
106 const scalar Re,
107 const scalar muc
108) const
109{
110 forceSuSp value(Zero);
111
112 const interpolation<vector>& HdotGradHInterp = *HdotGradHInterpPtr_;
113
114 value.Su()=
115 mass*3.0*constant::electromagnetic::mu0.value()/p.rho()
116 *magneticSusceptibility_/(magneticSusceptibility_ + 3)
117 *HdotGradHInterp.interpolate(p.coordinates(), p.currentTetIndices());
118
119 return value;
120}
121
122
123// ************************************************************************* //
Calculates particle paramagnetic (magnetic field) force.
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.
virtual ~ParamagneticForce()
Destructor.
ParamagneticForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Abstract base class for particle forces.
const fvMesh & mesh() const noexcept
Return the mesh database.
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
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition forceSuSpI.H:54
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)
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
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.
void deleteDemandDrivenData(DataPtr &dataPtr)
dictionary dict