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 | Copyright (C) 2011 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_(NULL),
46  magneticSusceptibility_
47  (
48  readScalar(this->coeffs().lookup("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 scalar dt,
101  const scalar mass,
102  const scalar Re,
103  const scalar muc
104 ) const
105 {
106  forceSuSp value(vector::zero, 0.0);
107 
108  const interpolation<vector>& HdotGradHInterp = *HdotGradHInterpPtr_;
109 
110  value.Su()=
111  mass*3.0*constant::electromagnetic::mu0.value()/p.rho()
112  *magneticSusceptibility_/(magneticSusceptibility_ + 3)
113  *HdotGradHInterp.interpolate(p.position(), p.currentTetIndices());
114 
115  return value;
116 }
117 
118 
119 // ************************************************************************* //
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
p
p
Definition: pEqn.H:62
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::ParamagneticForce::~ParamagneticForce
virtual ~ParamagneticForce()
Destructor.
Definition: ParamagneticForce.C:69
Foam::ParamagneticForce::HdotGradHInterpPtr_
const interpolation< vector > * HdotGradHInterpPtr_
HdotGradH interpolator.
Definition: ParamagneticForce.H:63
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
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.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
ParamagneticForce.H
Foam::ParamagneticForce::ParamagneticForce
ParamagneticForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Definition: ParamagneticForce.C:34
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::forceSuSp
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:58
Foam::ParticleForce
Abstract base class for particle forces.
Definition: ParticleForce.H:53
Foam::ParamagneticForce::calcNonCoupled
virtual forceSuSp calcNonCoupled(const typename CloudType::parcelType &p, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
Definition: ParamagneticForce.C:98
Foam::ParamagneticForce
Calculates particle paramagnetic (magnetic field) force.
Definition: ParamagneticForce.H:53
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::ParamagneticForce::cacheFields
virtual void cacheFields(const bool store)
Cache fields.
Definition: ParamagneticForce.C:76
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
readScalar
#define readScalar
Definition: doubleScalar.C:38
electromagneticConstants.H
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:217
Foam::ParamagneticForce::HdotGradHName_
const word HdotGradHName_
Name of paramagnetic field strength field - default = "HdotGradH".
Definition: ParamagneticForce.H:60
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::forceSuSp::Su
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:56
Foam::ParamagneticForce::magneticSusceptibility_
const scalar magneticSusceptibility_
Magnetic susceptibility of particle.
Definition: ParamagneticForce.H:66
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress