ParamagneticForce.H
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-------------------------------------------------------------------------------
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
26Class
27 Foam::ParamagneticForce
28
29Group
30 grpLagrangianIntermediateForceSubModels
31
32Description
33 Calculates particle paramagnetic (magnetic field) force
34
35SourceFiles
36 ParamagneticForceI.H
37 ParamagneticForce.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef ParamagneticForce_H
42#define ParamagneticForce_H
43
44#include "ParticleForce.H"
45#include "interpolation.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52class fvMesh;
53
54/*---------------------------------------------------------------------------*\
55 Class ParamagneticForce Declaration
56\*---------------------------------------------------------------------------*/
57
58template<class CloudType>
60:
61 public ParticleForce<CloudType>
62{
63 // Private data
64
65 //- Name of paramagnetic field strength field - default = "HdotGradH"
66 const word HdotGradHName_;
67
68 //- HdotGradH interpolator
69 const interpolation<vector>* HdotGradHInterpPtr_;
70
71 //- Magnetic susceptibility of particle
72 const scalar magneticSusceptibility_;
73
74
75public:
76
77 //- Runtime type information
78 TypeName("paramagnetic");
79
80
81 // Constructors
82
83 //- Construct from mesh
85 (
87 const fvMesh& mesh,
88 const dictionary& dict
89 );
90
91 //- Construct copy
93
94 //- Construct and return a clone
96 {
98 (
100 );
101 }
102
103
104 //- Destructor
105 virtual ~ParamagneticForce();
106
107
108 // Member Functions
109
110 // Access
111
112 //- Return the name of paramagnetic field strength field
113 const word& HdotGradHName() const;
114
115 //- Return the magnetic susceptibility of particle
116 scalar magneticSusceptibility() const;
117
118
119 // Evaluation
120
121 //- Cache fields
122 virtual void cacheFields(const bool store);
123
124 //- Calculate the non-coupled force
126 (
127 const typename CloudType::parcelType& p,
128 const typename CloudType::parcelType::trackingData& td,
129 const scalar dt,
130 const scalar mass,
131 const scalar Re,
132 const scalar muc
133 ) const;
134};
135
136
137// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138
139} // End namespace Foam
140
141// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142
143#include "ParamagneticForceI.H"
144
145#ifdef NoRepository
146 #include "ParamagneticForce.C"
147#endif
148
149// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150
151#endif
152
153// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Calculates particle paramagnetic (magnetic field) force.
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.
const word & HdotGradHName() const
Return the name of paramagnetic field strength field.
scalar magneticSusceptibility() const
Return the magnetic susceptibility of particle.
TypeName("paramagnetic")
Runtime type information.
virtual ~ParamagneticForce()
Destructor.
Abstract base class for particle forces.
Definition: ParticleForce.H:60
const CloudType & owner() const
Return const access to the cloud owner.
const fvMesh & mesh() const
Return the mesh database.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:67
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
Namespace for OpenFOAM.
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73