ParamagneticForce.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "ParamagneticForce.H"
27 #include "demandDrivenData.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  CloudType& owner,
36  const fvMesh& mesh,
37  const dictionary& dict
38 )
39 :
40  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
41  HdotGradHName_
42  (
43  this->coeffs().template lookupOrDefault<word>("HdotGradH", "HdotGradH")
44  ),
45  HdotGradHInterpPtr_(nullptr),
46  magneticSusceptibility_
47  (
48  this->coeffs().template lookup<scalar>("magneticSusceptibility")
49  )
50 {}
51 
52 
53 template<class CloudType>
55 (
56  const ParamagneticForce& pf
57 )
58 :
60  HdotGradHName_(pf.HdotGradHName_),
61  HdotGradHInterpPtr_(pf.HdotGradHInterpPtr_),
62  magneticSusceptibility_(pf.magneticSusceptibility_)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
67 
68 template<class CloudType>
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class CloudType>
77 {
78  if (store)
79  {
80  const volVectorField& HdotGradH =
81  this->mesh().template lookupObject<volVectorField>(HdotGradHName_);
82 
83  HdotGradHInterpPtr_ = interpolation<vector>::New
84  (
85  this->owner().solution().interpolationSchemes(),
86  HdotGradH
87  ).ptr();
88  }
89  else
90  {
91  deleteDemandDrivenData(HdotGradHInterpPtr_);
92  }
93 }
94 
95 
96 template<class CloudType>
98 (
99  const typename CloudType::parcelType& p,
100  const typename CloudType::parcelType::trackingData& td,
101  const scalar dt,
102  const scalar mass,
103  const scalar Re,
104  const scalar muc
105 ) const
106 {
107  forceSuSp value(Zero, 0.0);
108 
109  const interpolation<vector>& HdotGradHInterp = *HdotGradHInterpPtr_;
110 
111  value.Su()=
112  mass*3.0*constant::electromagnetic::mu0.value()/p.rho()
113  *magneticSusceptibility_/(magneticSusceptibility_ + 3)
114  *HdotGradHInterp.interpolate
115  (
116  p.coordinates(),
117  p.currentTetIndices(td.mesh)
118  );
119 
120  return value;
121 }
122 
123 
124 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Generic GeometricField class.
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.
Definition: ParticleForce.H:54
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:64
const vector & Su() const
Return const access to the explicit contribution [kg m/s^2].
Definition: forceSuSpI.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract base class for interpolation.
Definition: interpolation.H:55
static autoPtr< interpolation< Type > > New(const word &interpolationType, const VolField< Type > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:56
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:51
A class for handling words, derived from string.
Definition: word.H:62
Template functions to aid in the implementation of demand driven data.
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
static const zero Zero
Definition: zero.H:97
void deleteDemandDrivenData(DataType *&dataPtr)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict
volScalarField & p