ParamagneticForce.H
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-2020 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 Class
25  Foam::ParamagneticForce
26 
27 Description
28  Calculates particle paramagnetic (magnetic field) force
29 
30 SourceFiles
31  ParamagneticForceI.H
32  ParamagneticForce.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef ParamagneticForce_H
37 #define ParamagneticForce_H
38 
39 #include "ParticleForce.H"
40 #include "interpolation.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 class fvMesh;
48 
49 /*---------------------------------------------------------------------------*\
50  Class ParamagneticForce Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class CloudType>
55 :
56  public ParticleForce<CloudType>
57 {
58  // Private Data
59 
60  //- Name of paramagnetic field strength field - default = "HdotGradH"
61  const word HdotGradHName_;
62 
63  //- HdotGradH interpolator
64  const interpolation<vector>* HdotGradHInterpPtr_;
65 
66  //- Magnetic susceptibility of particle
67  const scalar magneticSusceptibility_;
68 
69 
70 public:
71 
72  //- Runtime type information
73  TypeName("paramagnetic");
74 
75 
76  // Constructors
77 
78  //- Construct from mesh
80  (
82  const fvMesh& mesh,
83  const dictionary& dict
84  );
85 
86  //- Construct copy
88 
89  //- Construct and return a clone
91  {
93  (
95  );
96  }
97 
98 
99  //- Destructor
100  virtual ~ParamagneticForce();
101 
102 
103  // Member Functions
104 
105  // Access
106 
107  //- Return the name of paramagnetic field strength field
108  const word& HdotGradHName() const;
109 
110  //- Return the magnetic susceptibility of particle
111  scalar magneticSusceptibility() const;
112 
113 
114  // Evaluation
115 
116  //- Cache fields
117  virtual void cacheFields(const bool store);
118 
119  //- Calculate the non-coupled force
120  virtual forceSuSp calcNonCoupled
121  (
122  const typename CloudType::parcelType& p,
123  const typename CloudType::parcelType::trackingData& td,
124  const scalar dt,
125  const scalar mass,
126  const scalar Re,
127  const scalar muc
128  ) const;
129 };
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 } // End namespace Foam
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 #include "ParamagneticForceI.H"
139 
140 #ifdef NoRepository
141  #include "ParamagneticForce.C"
142 #endif
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 #endif
147 
148 // ************************************************************************* //
dictionary dict
virtual void cacheFields(const bool store)
Cache fields.
TypeName("paramagnetic")
Runtime type information.
Calculates particle paramagnetic (magnetic field) force.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
volScalarField muc(IOobject(IOobject::groupName("mu", continuousPhaseName), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), rhoc *continuousPhaseViscosity->nu())
scalar magneticSusceptibility() const
Return the magnetic susceptibility of particle.
ParamagneticForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
const fvMesh & mesh() const
Return the mesh database.
Abstract base class for particle forces.
Definition: ParticleForce.H:53
const CloudType & owner() const
Return const access to the cloud owner.
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.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:61
A class for handling words, derived from string.
Definition: word.H:59
virtual ~ParamagneticForce()
Destructor.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
const word & HdotGradHName() const
Return the name of paramagnetic field strength field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Abstract base class for interpolation.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
virtual autoPtr< ParticleForce< CloudType > > clone() const
Construct and return a clone.
Namespace for OpenFOAM.