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 -------------------------------------------------------------------------------
10 License
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 "ParamagneticForce.H"
29 #include "demandDrivenData.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  CloudType& owner,
38  const fvMesh& mesh,
39  const dictionary& dict
40 )
41 :
42  ParticleForce<CloudType>(owner, mesh, dict, typeName, true),
43  HdotGradHName_
44  (
45  this->coeffs().template lookupOrDefault<word>("HdotGradH", "HdotGradH")
46  ),
47  HdotGradHInterpPtr_(nullptr),
48  magneticSusceptibility_
49  (
50  this->coeffs().getScalar("magneticSusceptibility")
51  )
52 {}
53 
54 
55 template<class CloudType>
57 (
58  const ParamagneticForce& pf
59 )
60 :
62  HdotGradHName_(pf.HdotGradHName_),
63  HdotGradHInterpPtr_(pf.HdotGradHInterpPtr_),
64  magneticSusceptibility_(pf.magneticSusceptibility_)
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
69 
70 template<class CloudType>
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
77 template<class CloudType>
79 {
80  if (store)
81  {
82  const volVectorField& HdotGradH =
83  this->mesh().template lookupObject<volVectorField>(HdotGradHName_);
84 
85  HdotGradHInterpPtr_ = interpolation<vector>::New
86  (
87  this->owner().solution().interpolationSchemes(),
88  HdotGradH
89  ).ptr();
90  }
91  else
92  {
93  deleteDemandDrivenData(HdotGradHInterpPtr_);
94  }
95 }
96 
97 
98 template<class CloudType>
100 (
101  const typename CloudType::parcelType& p,
102  const typename CloudType::parcelType::trackingData& td,
103  const scalar dt,
104  const scalar mass,
105  const scalar Re,
106  const scalar muc
107 ) const
108 {
109  forceSuSp value(Zero);
110 
111  const interpolation<vector>& HdotGradHInterp = *HdotGradHInterpPtr_;
112 
113  value.Su()=
114  mass*3.0*constant::electromagnetic::mu0.value()/p.rho()
115  *magneticSusceptibility_/(magneticSusceptibility_ + 3)
116  *HdotGradHInterp.interpolate(p.coordinates(), p.currentTetIndices());
117 
118  return value;
119 }
120 
121 
122 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::ParamagneticForce::calcNonCoupled
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.
Definition: ParamagneticForce.C:100
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:50
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::ParamagneticForce::~ParamagneticForce
virtual ~ParamagneticForce()
Destructor.
Definition: ParamagneticForce.C:71
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::forceSuSp::Su
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:64
ParamagneticForce.H
Foam::ParamagneticForce::ParamagneticForce
ParamagneticForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Definition: ParamagneticForce.C:36
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::forceSuSp
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:63
Foam::ParticleForce
Abstract base class for particle forces.
Definition: ParticleForce.H:59
Foam::ParamagneticForce
Calculates particle paramagnetic (magnetic field) force.
Definition: ParamagneticForce.H:58
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:95
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::ParamagneticForce::cacheFields
virtual void cacheFields(const bool store)
Cache fields.
Definition: ParamagneticForce.C:78
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
electromagneticConstants.H
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::constant::electromagnetic::mu0
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::interpolation::interpolate
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.