Loading...
Searching...
No Matches
moleculeI.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30
32:
33 siteReferencePositions_(),
34 siteMasses_(),
35 siteCharges_(),
36 siteIds_(),
37 pairPotentialSites_(),
38 electrostaticSites_(),
39 momentOfInertia_(Zero),
40 mass_(0)
41{}
42
43
45(
46 const dictionary& dict
47)
48:
49 siteReferencePositions_(dict.lookup("siteReferencePositions")),
50 siteMasses_(dict.lookup("siteMasses")),
51 siteCharges_(dict.lookup("siteCharges")),
52 siteIds_(List<word>(dict.lookup("siteIds")).size(), -1),
53 pairPotentialSites_(),
54 electrostaticSites_(),
55 momentOfInertia_(),
56 mass_()
57{
58 checkSiteListSizes();
59
60 setInteracionSiteBools
61 (
62 List<word>(dict.lookup("siteIds")),
63 List<word>(dict.lookup("pairPotentialSiteIds"))
64 );
65
66 mass_ = sum(siteMasses_);
67
68 vector centreOfMass(Zero);
69
70 // Calculate the centre of mass of the body and subtract it from each
71 // position
72
73 forAll(siteReferencePositions_, i)
74 {
75 centreOfMass += siteReferencePositions_[i]*siteMasses_[i];
76 }
77
78 centreOfMass /= mass_;
79
80 siteReferencePositions_ -= centreOfMass;
81
82 if (siteIds_.size() == 1)
83 {
84 // Single site molecule - no rotational motion.
85
86 siteReferencePositions_[0] = Zero;
87
88 momentOfInertia_ = diagTensor(-1, -1, -1);
89 }
90 else if (linearMoleculeTest())
91 {
92 // Linear molecule.
93
94 Info<< nl << "Linear molecule." << endl;
95
96 const vector dir =
98 (
99 siteReferencePositions_[1] - siteReferencePositions_[0]
100 );
101
102 tensor Q = rotationTensor(dir, vector(1,0,0));
103
104 siteReferencePositions_ = (Q & siteReferencePositions_);
105
106 // The rotation was around the centre of mass but remove any
107 // components that have crept in due to floating point errors
108
109 centreOfMass = Zero;
110
111 forAll(siteReferencePositions_, i)
112 {
113 centreOfMass += siteReferencePositions_[i]*siteMasses_[i];
114 }
115
116 centreOfMass /= mass_;
117
118 siteReferencePositions_ -= centreOfMass;
119
120 diagTensor momOfInertia = Zero;
121
122 forAll(siteReferencePositions_, i)
123 {
124 const vector& p(siteReferencePositions_[i]);
125
126 momOfInertia +=
127 siteMasses_[i]*diagTensor(0, p.x()*p.x(), p.x()*p.x());
128 }
129
130 momentOfInertia_ = diagTensor
131 (
132 -1,
133 momOfInertia.yy(),
134 momOfInertia.zz()
135 );
136 }
137 else
138 {
139 // Fully 6DOF molecule
140
141 // Calculate the inertia tensor in the current orientation
142
143 symmTensor momOfInertia(Zero);
144
145 forAll(siteReferencePositions_, i)
146 {
147 const vector& p(siteReferencePositions_[i]);
148
149 momOfInertia += siteMasses_[i]*symmTensor
150 (
151 p.y()*p.y() + p.z()*p.z(), -p.x()*p.y(), -p.x()*p.z(),
152 p.x()*p.x() + p.z()*p.z(), -p.y()*p.z(),
153 p.x()*p.x() + p.y()*p.y()
154 );
155 }
156
157 if (eigenValues(momOfInertia).x() < VSMALL)
158 {
160 << "An eigenvalue of the inertia tensor is zero, the molecule "
161 << dict.name() << " is not a valid 6DOF shape."
162 << nl << abort(FatalError);
163 }
164
165 // Normalise the inertia tensor magnitude to avoid SMALL numbers in the
166 // components causing problems
167
168 momOfInertia /= eigenValues(momOfInertia).x();
169
170 tensor e(eigenVectors(momOfInertia));
171
172 // Calculate the transformation between the principle axes to the
173 // global axes
174
175 tensor Q =
176 vector(1,0,0)*e.x() + vector(0,1,0)*e.y() + vector(0,0,1)*e.z();
177
178 // Transform the site positions
179
180 siteReferencePositions_ = (Q & siteReferencePositions_);
181
182 // Recalculating the moment of inertia with the new site positions
183
184 // The rotation was around the centre of mass but remove any
185 // components that have crept in due to floating point errors
186
187 centreOfMass = Zero;
188
189 forAll(siteReferencePositions_, i)
190 {
191 centreOfMass += siteReferencePositions_[i]*siteMasses_[i];
192 }
193
194 centreOfMass /= mass_;
195
196 siteReferencePositions_ -= centreOfMass;
197
198 // Calculate the moment of inertia in the principle component
199 // reference frame
200
201 momOfInertia = Zero;
202
203 forAll(siteReferencePositions_, i)
204 {
205 const vector& p(siteReferencePositions_[i]);
206
207 momOfInertia += siteMasses_[i]*symmTensor
208 (
209 p.y()*p.y() + p.z()*p.z(), -p.x()*p.y(), -p.x()*p.z(),
210 p.x()*p.x() + p.z()*p.z(), -p.y()*p.z(),
211 p.x()*p.x() + p.y()*p.y()
212 );
213 }
214
215 momentOfInertia_ = diagTensor
216 (
217 momOfInertia.xx(),
218 momOfInertia.yy(),
219 momOfInertia.zz()
220 );
221 }
222}
223
224
226(
227 const polyMesh& mesh,
229 const label celli,
230 const label tetFacei,
231 const label tetPti,
232 const tensor& Q,
233 const vector& v,
234 const vector& a,
235 const vector& pi,
236 const vector& tau,
237 const vector& specialPosition,
238 const constantProperties& constProps,
239 const label special,
240 const label id
241
242)
243:
244 particle(mesh, coordinates, celli, tetFacei, tetPti),
245 Q_(Q),
246 v_(v),
247 a_(a),
248 pi_(pi),
249 tau_(tau),
250 specialPosition_(specialPosition),
251 potentialEnergy_(0.0),
252 rf_(Zero),
253 special_(special),
254 id_(id),
255 siteForces_(constProps.nSites(), Zero),
256 sitePositions_(constProps.nSites())
257{
258 setSitePositions(constProps);
259}
260
261
263(
264 const polyMesh& mesh,
265 const vector& position,
266 const label celli,
267 const tensor& Q,
268 const vector& v,
269 const vector& a,
270 const vector& pi,
271 const vector& tau,
272 const vector& specialPosition,
273 const constantProperties& constProps,
274 const label special,
275 const label id
276
277)
278:
279 particle(mesh, position, celli),
280 Q_(Q),
281 v_(v),
282 a_(a),
283 pi_(pi),
284 tau_(tau),
285 specialPosition_(specialPosition),
286 potentialEnergy_(0.0),
287 rf_(Zero),
288 special_(special),
289 id_(id),
290 siteForces_(constProps.nSites(), Zero),
291 sitePositions_(constProps.nSites())
292{
293 setSitePositions(constProps);
294}
295
296
297// * * * * * * * constantProperties Private Member Functions * * * * * * * * //
298
299inline void Foam::molecule::constantProperties::checkSiteListSizes() const
300{
301 if
302 (
303 siteIds_.size() != siteReferencePositions_.size()
304 || siteIds_.size() != siteCharges_.size()
305 )
306 {
308 << "Sizes of site id, charge and "
309 << "referencePositions are not the same. "
310 << nl << abort(FatalError);
311 }
312}
313
314
315inline void Foam::molecule::constantProperties::setInteracionSiteBools
316(
317 const List<word>& siteIds,
318 const List<word>& pairPotSiteIds
319)
320{
321 pairPotentialSites_.setSize(siteIds_.size());
322
323 electrostaticSites_.setSize(siteIds_.size());
324
325 forAll(siteIds_, i)
326 {
327 const word& id = siteIds[i];
328
329 pairPotentialSites_[i] = pairPotSiteIds.found(id);
330 electrostaticSites_[i] = (mag(siteCharges_[i]) > VSMALL);
331 }
332}
333
334
335inline bool Foam::molecule::constantProperties::linearMoleculeTest() const
336{
337 if (siteIds_.size() == 2)
338 {
339 return true;
340 }
341
342 const vector refDir =
344 (
345 siteReferencePositions_[1] - siteReferencePositions_[0]
346 );
347
348 for
349 (
350 label i = 2;
351 i < siteReferencePositions_.size();
352 i++
353 )
354 {
355 const vector dir =
357 (
358 siteReferencePositions_[i] - siteReferencePositions_[i-1]
359 );
360
361 if (mag(refDir & dir) < 1 - SMALL)
362 {
363 return false;
364 }
365 }
366
367 return true;
368}
369
370
371// * * * * * * * constantProperties Member Functions * * * * * * * * * * * * //
372
373inline const Foam::Field<Foam::vector>&
375{
376 return siteReferencePositions_;
377}
378
379
380inline const Foam::List<Foam::scalar>&
382{
383 return siteMasses_;
384}
385
386
387inline const Foam::List<Foam::scalar>&
389{
390 return siteCharges_;
391}
392
393
394inline const Foam::List<Foam::label>&
396{
397 return siteIds_;
398}
399
400
403{
404 return siteIds_;
405}
406
407
408inline const Foam::List<bool>&
410{
411 return pairPotentialSites_;
412}
413
414
416(
417 label sId
418) const
419{
420 label s = siteIds_.find(sId);
421
422 if (s == -1)
423 {
425 << sId << " site not found."
426 << nl << abort(FatalError);
428
429 return pairPotentialSites_[s];
430}
431
432
433inline const Foam::List<bool>&
435{
436 return electrostaticSites_;
437}
438
439
441(
442 label sId
443) const
444{
445 label s = siteIds_.find(sId);
446
447 if (s == -1)
448 {
450 << sId << " site not found."
451 << nl << abort(FatalError);
453
454 return electrostaticSites_[s];
455}
456
457
458inline const Foam::diagTensor&
460{
461 return momentOfInertia_;
462}
463
466{
467 return ((momentOfInertia_.xx() < 0) && (momentOfInertia_.yy() > 0));
468}
469
472{
473 return (momentOfInertia_.zz() < 0);
474}
475
476
478{
479 if (linearMolecule())
480 {
481 return 5;
482 }
483 else if (pointMolecule())
484 {
485 return 3;
486 }
487 else
488 {
489 return 6;
490 }
491}
492
494inline Foam::scalar Foam::molecule::constantProperties::mass() const
495{
496 return mass_;
497}
498
499
500inline Foam::label Foam::molecule::constantProperties::nSites() const
501{
502 return siteIds_.size();
503
504}
505
506
507// * * * * * * * * * * * * molecule Member Functions * * * * * * * * * * * * //
509inline const Foam::tensor& Foam::molecule::Q() const
510{
511 return Q_;
512}
513
516{
517 return Q_;
518}
519
521inline const Foam::vector& Foam::molecule::v() const
522{
523 return v_;
524}
525
528{
529 return v_;
530}
531
533inline const Foam::vector& Foam::molecule::a() const
534{
535 return a_;
536}
537
540{
541 return a_;
542}
543
545inline const Foam::vector& Foam::molecule::pi() const
546{
547 return pi_;
548}
549
552{
553 return pi_;
554}
555
557inline const Foam::vector& Foam::molecule::tau() const
558{
559 return tau_;
560}
561
564{
565 return tau_;
566}
567
570{
571 return siteForces_;
572}
573
576{
577 return siteForces_;
578}
579
582{
583 return sitePositions_;
584}
585
588{
589 return sitePositions_;
590}
591
594{
595 return specialPosition_;
596}
597
600{
601 return specialPosition_;
602}
603
605inline Foam::scalar Foam::molecule::potentialEnergy() const
606{
607 return potentialEnergy_;
608}
609
611inline Foam::scalar& Foam::molecule::potentialEnergy()
612{
613 return potentialEnergy_;
614}
615
617inline const Foam::tensor& Foam::molecule::rf() const
618{
619 return rf_;
620}
621
624{
625 return rf_;
626}
627
629inline Foam::label Foam::molecule::special() const
630{
631 return special_;
632}
633
635inline bool Foam::molecule::tethered() const
636{
637 return special_ == SPECIAL_TETHERED;
638}
639
640
641inline Foam::label Foam::molecule::id() const
642{
643 return id_;
644}
645
646
647// ************************************************************************* //
constexpr scalar pi(M_PI)
const Cmpt & yy() const noexcept
Definition DiagTensor.H:133
const Cmpt & zz() const noexcept
Definition DiagTensor.H:134
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
const Cmpt & yy() const noexcept
Definition SymmTensor.H:154
const Cmpt & xx() const noexcept
Definition SymmTensor.H:150
const Cmpt & zz() const noexcept
Definition SymmTensor.H:158
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Class to hold molecule constant properties.
Definition molecule.H:91
const diagTensor & momentOfInertia() const
Definition moleculeI.H:452
const Field< vector > & siteReferencePositions() const
Definition moleculeI.H:367
const List< label > & siteIds() const
Definition moleculeI.H:388
const List< bool > & pairPotentialSites() const
Definition moleculeI.H:402
const List< scalar > & siteCharges() const
Definition moleculeI.H:381
bool pairPotentialSite(label sId) const
Definition moleculeI.H:409
bool electrostaticSite(label sId) const
Definition moleculeI.H:434
const List< scalar > & siteMasses() const
Definition moleculeI.H:374
const List< bool > & electrostaticSites() const
Definition moleculeI.H:427
const tensor & Q() const
Definition moleculeI.H:502
const vector & v() const
Definition moleculeI.H:514
const List< vector > & sitePositions() const
Definition moleculeI.H:574
const vector & a() const
Definition moleculeI.H:526
const List< vector > & siteForces() const
Definition moleculeI.H:562
const vector & specialPosition() const
Definition moleculeI.H:586
void setSitePositions(const constantProperties &constProps)
Definition molecule.C:219
label special() const
Definition moleculeI.H:622
const tensor & rf() const
Definition moleculeI.H:610
const vector & tau() const
Definition moleculeI.H:550
molecule(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const tensor &Q, const vector &v, const vector &a, const vector &pi, const vector &tau, const vector &specialPosition, const constantProperties &constProps, const label special, const label id)
Construct from components.
Definition moleculeI.H:219
scalar potentialEnergy() const
Definition moleculeI.H:598
const vector & pi() const
Definition moleculeI.H:538
bool tethered() const
Definition moleculeI.H:628
label id() const
Definition moleculeI.H:634
Base particle class.
Definition particle.H:72
vector position() const
Return current particle position.
Definition particleI.H:283
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Lookup type of boundary radiation properties.
Definition lookup.H:60
Tensor of scalars, i.e. Tensor<scalar>.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
DiagTensor< scalar > diagTensor
DiagTensor of scalars, i.e. DiagTensor<scalar>.
Definition diagTensor.H:49
dimensionedTensor eigenVectors(const dimensionedSymmTensor &dt)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
messageStream Info
Information stream (stdout output on master, null elsewhere).
Tensor< scalar > tensor
Definition symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299