moleculeI.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-2024 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 "molecule.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 :
32  siteReferencePositions_(Field<vector>(0)),
33  siteMasses_(List<scalar>(0)),
34  siteCharges_(List<scalar>(0)),
35  siteIds_(List<label>(0)),
36  pairPotentialSites_(List<bool>(false)),
37  electrostaticSites_(List<bool>(false)),
38  momentOfInertia_(diagTensor(0, 0, 0)),
39  mass_(0)
40 {}
41 
42 
44 (
45  const dictionary& dict
46 )
47 :
48  siteReferencePositions_(dict.lookup("siteReferencePositions")),
49  siteMasses_(dict.lookup("siteMasses")),
50  siteCharges_(dict.lookup("siteCharges")),
51  siteIds_(List<word>(dict.lookup("siteIds")).size(), -1),
52  pairPotentialSites_(),
53  electrostaticSites_(),
54  momentOfInertia_(),
55  mass_()
56 {
57  checkSiteListSizes();
58 
59  setInteractionSiteBools
60  (
61  List<word>(dict.lookup("siteIds")),
62  List<word>(dict.lookup("pairPotentialSiteIds"))
63  );
64 
65  mass_ = sum(siteMasses_);
66 
67  vector centreOfMass(Zero);
68 
69  // Calculate the centre of mass of the body and subtract it from each
70  // position
71 
72  forAll(siteReferencePositions_, i)
73  {
74  centreOfMass += siteReferencePositions_[i]*siteMasses_[i];
75  }
76 
77  centreOfMass /= mass_;
78 
79  siteReferencePositions_ -= centreOfMass;
80 
81  if (siteIds_.size() == 1)
82  {
83  // Single site molecule - no rotational motion.
84 
85  siteReferencePositions_[0] = Zero;
86 
87  momentOfInertia_ = diagTensor(-1, -1, -1);
88  }
89  else if (linearMoleculeTest())
90  {
91  // Linear molecule.
92 
93  Info<< nl << "Linear molecule." << endl;
94 
95  vector dir = siteReferencePositions_[1] - siteReferencePositions_[0];
96 
97  dir /= mag(dir);
98 
99  tensor Q = rotationTensor(dir, vector(1,0,0));
100 
101  siteReferencePositions_ = (Q & siteReferencePositions_);
102 
103  // The rotation was around the centre of mass but remove any
104  // components that have crept in due to floating point errors
105 
106  centreOfMass = Zero;
107 
108  forAll(siteReferencePositions_, i)
109  {
110  centreOfMass += siteReferencePositions_[i]*siteMasses_[i];
111  }
112 
113  centreOfMass /= mass_;
114 
115  siteReferencePositions_ -= centreOfMass;
116 
117  diagTensor momOfInertia = Zero;
118 
119  forAll(siteReferencePositions_, i)
120  {
121  const vector& p(siteReferencePositions_[i]);
122 
123  momOfInertia +=
124  siteMasses_[i]*diagTensor(0, p.x()*p.x(), p.x()*p.x());
125  }
126 
127  momentOfInertia_ = diagTensor
128  (
129  -1,
130  momOfInertia.yy(),
131  momOfInertia.zz()
132  );
133  }
134  else
135  {
136  // Fully 6DOF molecule
137 
138  // Calculate the inertia tensor in the current orientation
139 
140  tensor momOfInertia(Zero);
141 
142  forAll(siteReferencePositions_, i)
143  {
144  const vector& p(siteReferencePositions_[i]);
145 
146  momOfInertia += siteMasses_[i]*tensor
147  (
148  p.y()*p.y() + p.z()*p.z(), -p.x()*p.y(), -p.x()*p.z(),
149  -p.y()*p.x(), p.x()*p.x() + p.z()*p.z(), -p.y()*p.z(),
150  -p.z()*p.x(), -p.z()*p.y(), p.x()*p.x() + p.y()*p.y()
151  );
152  }
153 
154  if (eigenValues(momOfInertia).x() < vSmall)
155  {
157  << "An eigenvalue of the inertia tensor is zero, the molecule "
158  << dict.name() << " is not a valid 6DOF shape."
159  << nl << abort(FatalError);
160  }
161 
162  // Normalise the inertia tensor magnitude to avoid small numbers in the
163  // components causing problems
164 
165  momOfInertia /= eigenValues(momOfInertia).x();
166 
167  tensor e = eigenVectors(momOfInertia);
168 
169  // Calculate the transformation between the principle axes to the
170  // global axes
171 
172  tensor Q =
173  vector(1,0,0)*e.x() + vector(0,1,0)*e.y() + vector(0,0,1)*e.z();
174 
175  // Transform the site positions
176 
177  siteReferencePositions_ = (Q & siteReferencePositions_);
178 
179  // Recalculating the moment of inertia with the new site positions
180 
181  // The rotation was around the centre of mass but remove any
182  // components that have crept in due to floating point errors
183 
184  centreOfMass = Zero;
185 
186  forAll(siteReferencePositions_, i)
187  {
188  centreOfMass += siteReferencePositions_[i]*siteMasses_[i];
189  }
190 
191  centreOfMass /= mass_;
192 
193  siteReferencePositions_ -= centreOfMass;
194 
195  // Calculate the moment of inertia in the principle component
196  // reference frame
197 
198  momOfInertia = Zero;
199 
200  forAll(siteReferencePositions_, i)
201  {
202  const vector& p(siteReferencePositions_[i]);
203 
204  momOfInertia += siteMasses_[i]*tensor
205  (
206  p.y()*p.y() + p.z()*p.z(), -p.x()*p.y(), -p.x()*p.z(),
207  -p.y()*p.x(), p.x()*p.x() + p.z()*p.z(), -p.y()*p.z(),
208  -p.z()*p.x(), -p.z()*p.y(), p.x()*p.x() + p.y()*p.y()
209  );
210  }
211 
212  momentOfInertia_ = diagTensor
213  (
214  momOfInertia.xx(),
215  momOfInertia.yy(),
216  momOfInertia.zz()
217  );
218  }
219 }
220 
221 
223 (
224  const polyMesh& mesh,
225  const vector& position,
226  const label celli,
227  label& nLocateBoundaryHits,
228  const tensor& Q,
229  const vector& v,
230  const vector& a,
231  const vector& pi,
232  const vector& tau,
233  const vector& specialPosition,
234  const constantProperties& constProps,
235  const label special,
236  const label id
237 
238 )
239 :
240  particle(mesh, position, celli, nLocateBoundaryHits),
241  Q_(Q),
242  v_(v),
243  a_(a),
244  pi_(pi),
245  tau_(tau),
246  specialPosition_(specialPosition),
247  potentialEnergy_(0.0),
248  rf_(Zero),
249  special_(special),
250  id_(id),
251  siteForces_(constProps.nSites(), Zero),
252  sitePositions_(constProps.nSites())
253 {
254  setSitePositions(mesh, constProps);
255 }
256 
257 
258 // * * * * * * * constantProperties Private Member Functions * * * * * * * * //
259 
260 inline void Foam::molecule::constantProperties::checkSiteListSizes() const
261 {
262  if
263  (
264  siteIds_.size() != siteReferencePositions_.size()
265  || siteIds_.size() != siteCharges_.size()
266  )
267  {
269  << "Sizes of site id, charge and "
270  << "referencePositions are not the same. "
271  << nl << abort(FatalError);
272  }
273 }
274 
275 
276 inline void Foam::molecule::constantProperties::setInteractionSiteBools
277 (
278  const List<word>& siteIds,
279  const List<word>& pairPotSiteIds
280 )
281 {
282  pairPotentialSites_.setSize(siteIds_.size());
283 
284  electrostaticSites_.setSize(siteIds_.size());
285 
286  forAll(siteIds_, i)
287  {
288  const word& id(siteIds[i]);
289 
290  pairPotentialSites_[i] = (findIndex(pairPotSiteIds, id) > -1);
291 
292  electrostaticSites_[i] = (mag(siteCharges_[i]) > vSmall);
293  }
294 }
295 
296 
297 inline bool Foam::molecule::constantProperties::linearMoleculeTest() const
298 {
299  if (siteIds_.size() == 2)
300  {
301  return true;
302  }
303 
304  vector refDir = siteReferencePositions_[1] - siteReferencePositions_[0];
305 
306  refDir /= mag(refDir);
307 
308  for
309  (
310  label i = 2;
311  i < siteReferencePositions_.size();
312  i++
313  )
314  {
315  vector dir = siteReferencePositions_[i] - siteReferencePositions_[i-1];
316 
317  dir /= mag(dir);
318 
319  if (mag(refDir & dir) < 1 - small)
320  {
321  return false;
322  }
323  }
324 
325  return true;
326 }
327 
328 
329 // * * * * * * * constantProperties Member Functions * * * * * * * * * * * * //
330 
331 inline const Foam::Field<Foam::vector>&
333 {
334  return siteReferencePositions_;
335 }
336 
337 
338 inline const Foam::List<Foam::scalar>&
340 {
341  return siteMasses_;
342 }
343 
344 
345 inline const Foam::List<Foam::scalar>&
347 {
348  return siteCharges_;
349 }
350 
351 
352 inline const Foam::List<Foam::label>&
354 {
355  return siteIds_;
356 }
357 
358 
361 {
362  return siteIds_;
363 }
364 
365 
366 inline const Foam::List<bool>&
368 {
369  return pairPotentialSites_;
370 }
371 
372 
374 (
375  label sId
376 ) const
377 {
378  label s = findIndex(siteIds_, sId);
379 
380  if (s == -1)
381  {
383  << sId << " site not found."
384  << nl << abort(FatalError);
385  }
386 
387  return pairPotentialSites_[s];
388 }
389 
390 
391 inline const Foam::List<bool>&
393 {
394  return electrostaticSites_;
395 }
396 
397 
399 (
400  label sId
401 ) const
402 {
403  label s = findIndex(siteIds_, sId);
404 
405  if (s == -1)
406  {
408  << sId << " site not found."
409  << nl << abort(FatalError);
410  }
411 
412  return electrostaticSites_[s];
413 }
414 
415 
416 inline const Foam::diagTensor&
418 {
419  return momentOfInertia_;
420 }
421 
422 
424 {
425  return ((momentOfInertia_.xx() < 0) && (momentOfInertia_.yy() > 0));
426 }
427 
428 
430 {
431  return (momentOfInertia_.zz() < 0);
432 }
433 
434 
436 {
437  if (linearMolecule())
438  {
439  return 5;
440  }
441  else if (pointMolecule())
442  {
443  return 3;
444  }
445  else
446  {
447  return 6;
448  }
449 }
450 
451 
452 inline Foam::scalar Foam::molecule::constantProperties::mass() const
453 {
454  return mass_;
455 }
456 
457 
459 {
460  return siteIds_.size();
461 
462 }
463 
464 
465 // * * * * * * * * * * * * molecule Member Functions * * * * * * * * * * * * //
466 
467 inline const Foam::tensor& Foam::molecule::Q() const
468 {
469  return Q_;
470 }
471 
472 
474 {
475  return Q_;
476 }
477 
478 
479 inline const Foam::vector& Foam::molecule::v() const
480 {
481  return v_;
482 }
483 
484 
486 {
487  return v_;
488 }
489 
490 
491 inline const Foam::vector& Foam::molecule::a() const
492 {
493  return a_;
494 }
495 
496 
498 {
499  return a_;
500 }
501 
502 
503 inline const Foam::vector& Foam::molecule::pi() const
504 {
505  return pi_;
506 }
507 
508 
510 {
511  return pi_;
512 }
513 
514 
515 inline const Foam::vector& Foam::molecule::tau() const
516 {
517  return tau_;
518 }
519 
520 
522 {
523  return tau_;
524 }
525 
526 
528 {
529  return siteForces_;
530 }
531 
532 
534 {
535  return siteForces_;
536 }
537 
538 
540 {
541  return sitePositions_;
542 }
543 
544 
546 {
547  return sitePositions_;
548 }
549 
550 
552 {
553  return specialPosition_;
554 }
555 
556 
558 {
559  return specialPosition_;
560 }
561 
562 
563 inline Foam::scalar Foam::molecule::potentialEnergy() const
564 {
565  return potentialEnergy_;
566 }
567 
568 
569 inline Foam::scalar& Foam::molecule::potentialEnergy()
570 {
571  return potentialEnergy_;
572 }
573 
574 
575 inline const Foam::tensor& Foam::molecule::rf() const
576 {
577  return rf_;
578 }
579 
580 
582 {
583  return rf_;
584 }
585 
586 
588 {
589  return special_;
590 }
591 
592 
593 inline bool Foam::molecule::tethered() const
594 {
595  return special_ == SPECIAL_TETHERED;
596 }
597 
598 
600 {
601  return id_;
602 }
603 
604 
605 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const Cmpt & zz() const
Definition: DiagTensorI.H:97
const Cmpt & yy() const
Definition: DiagTensorI.H:91
Pre-declare SubField and related Field type.
Definition: Field.H:83
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const Cmpt & xx() const
Definition: TensorI.H:163
const Cmpt & zz() const
Definition: TensorI.H:219
const Cmpt & yy() const
Definition: TensorI.H:191
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Class to hold molecule constant properties.
Definition: molecule.H:91
const diagTensor & momentOfInertia() const
Definition: moleculeI.H:417
const Field< vector > & siteReferencePositions() const
Definition: moleculeI.H:332
const List< label > & siteIds() const
Definition: moleculeI.H:353
const List< bool > & pairPotentialSites() const
Definition: moleculeI.H:367
const List< scalar > & siteCharges() const
Definition: moleculeI.H:346
bool pairPotentialSite(label sId) const
Definition: moleculeI.H:374
bool electrostaticSite(label sId) const
Definition: moleculeI.H:399
const List< scalar > & siteMasses() const
Definition: moleculeI.H:339
const List< bool > & electrostaticSites() const
Definition: moleculeI.H:392
const tensor & Q() const
Definition: moleculeI.H:467
molecule(const polyMesh &mesh, const vector &position, const label celli, label &nLocateBoundaryHits, 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 a position and a cell, searching for the rest of the.
Definition: moleculeI.H:223
const vector & v() const
Definition: moleculeI.H:479
const List< vector > & sitePositions() const
Definition: moleculeI.H:539
const vector & a() const
Definition: moleculeI.H:491
const List< vector > & siteForces() const
Definition: moleculeI.H:527
const vector & specialPosition() const
Definition: moleculeI.H:551
void setSitePositions(const polyMesh &mesh, const constantProperties &constProps)
Definition: molecule.C:218
label special() const
Definition: moleculeI.H:587
const tensor & rf() const
Definition: moleculeI.H:575
const vector & tau() const
Definition: moleculeI.H:515
scalar potentialEnergy() const
Definition: moleculeI.H:563
const vector & pi() const
Definition: moleculeI.H:503
bool tethered() const
Definition: moleculeI.H:593
label id() const
Definition: moleculeI.H:599
Base particle class.
Definition: particle.H:83
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:255
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from unit vector n1 to n2.
Definition: transform.H:47
dimensionedTensor eigenVectors(const dimensionedTensor &dt)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
DiagTensor< scalar > diagTensor
A scalar version of the templated DiagTensor.
Definition: diagTensor.H:48
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedVector eigenValues(const dimensionedTensor &dt)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p