Chombo + EB + MF  3.2
RealVect.H
Go to the documentation of this file.
1 #ifdef CH_LANG_CC
2 /*
3  * _______ __
4  * / ___/ / ___ __ _ / / ___
5  * / /__/ _ \/ _ \/ V \/ _ \/ _ \
6  * \___/_//_/\___/_/_/_/_.__/\___/
7  * Please refer to Copyright.txt, in Chombo's root directory.
8  */
9 #endif
10 
11 #ifndef _REALVECT_H_
12 #define _REALVECT_H_
13 
14 #include <cstddef>
15 #include <cstdlib>
16 #include <cstring>
17 #include <iostream>
18 
19 #include "SPACE.H"
20 #include "Misc.H"
21 #include "REAL.H"
22 #include "SPACE.H"
23 #include "IntVect.H"
24 #include "Vector.H"
25 #include "SPMD.H"
26 #include "IndexTM.H"
27 #include "NamespaceHeader.H"
28 
29 //template<typename T, int n>
30 //class IndexTM;
31 
32 /// A Real vector in SpaceDim-dimensional space
33 /**
34  The class RealVect is an implementation of a Real vector in a
35  SpaceDim-dimensional space.
36  RealVect values are accessed using the operator[] function, as for a normal
37  C++ array. In addition, the basic arithmetic operators have been overloaded
38  to implement scaling and translation operations.
39 */
40 
41 class RealVect
42 {
43 public:
44 
45  /**
46  \name Constructors and Accessors
47  */
48  /*@{*/
49 
50  ///
51  /**
52  Construct a RealVect whose components are zero.
53  */
54  RealVect ();
55 
56  explicit RealVect (const Vector<Real>&);
57 
58  ///
59  /**
60  Construct a RealVect given the specific values for its
61  coordinates. D_DECL6 is a macro that sets the constructor to
62  take CH_SPACEDIM arguments.
63  */
64  RealVect (D_DECL6(Real i, Real j, Real k,
65  Real l, Real m, Real n));
66 
67  ///
68  /**
69  The copy constructor.
70  */
71  RealVect (const RealVect& rhs);
72 
73  ///
74  /**
75  Construct a RealVect from an IndexTM<Real,SpaceDim>.
76  */
78 
79  ///
80  /**
81  Construct a RealVect from an IntVect by coercing each component
82  from <tt>int</tt> to Real.
83  */
84  RealVect (const IntVect & iv)
85  {
86  for (int d=0 ; d<SpaceDim ; ++d)
87  {
88  vect[d] = (Real)iv[d];
89  }
90  }
91 
92  ///
93  /**
94  The assignment operator.
95  */
96  RealVect& operator= (const RealVect& rhs);
97 
98  ///
99  /**
100  Returns a modifiable lvalue reference to the <i>i</i>'th coordinate of the
101  RealVect.
102  */
103  inline
104  Real& operator[] (int i);
105 
106  ///
107  /**
108  Returns the <i>i</i>'th coordinate of the RealVect.
109  */
110  inline
111  const Real& operator[] (int i) const;
112 
113  /*@}*/
114 
115  /**
116  \name Comparison Operators
117  */
118  /*@{*/
119 
120  ///
121  /**
122  Returns true if this RealVect is equivalent to argument RealVect. All
123  comparisons between analogous components must be satisfied.
124  */
125  bool operator== (const RealVect& p) const;
126 
127  ///
128  /**
129  Returns true if this RealVect is different from argument RealVect.
130  All comparisons between analogous components must be satisfied.
131  */
132  bool operator!= (const RealVect& p) const;
133 
134  ///
135  /**
136  Returns true if this RealVect is less than argument RealVect. All
137  comparisons between analogous components must be satisfied. Note
138  that, since the comparison is component-wise, it is possible for
139  an RealVect to be neither greater than, less than, nor equal to
140  another.
141  */
142  bool operator< (const RealVect& p) const;
143 
144  ///
145  /**
146  Returns true if this RealVect is less than or equal to argument
147  RealVect. All comparisons between analogous components must be
148  satisfied. Note that, since the comparison is component-wise, it
149  is possible for an RealVect to be neither greater than or equal
150  to, less than or equal to, nor equal to another.
151  */
152  bool operator<= (const RealVect& p) const;
153 
154  ///
155  /**
156  Returns true if this RealVect is greater than argument RealVect.
157  All comparisons between analogous components must be satisfied.
158  Note that, since the comparison is component-wise, it is possible
159  for an RealVect to be neither greater than, less than, nor equal
160  to another.
161  */
162  bool operator> (const RealVect& p) const;
163 
164  ///
165  /**
166  Returns true if this RealVect is greater than or equal to argument
167  RealVect. All comparisons between analogous components must be
168  satisfied. Note that, since the comparison is component-wise, it
169  is possible for an RealVect to be neither greater than or equal
170  to, less than or equal to, nor equal to another.
171  */
172 
173  bool operator>= (const RealVect& p) const;
174 
175  /*@}*/
176 
177  /**
178  \name Arithmetic Operators
179  */
180  /*@{*/
181 
182  ///
183  /**
184  Modifies this RealVect by addition of a scalar to each component.
185  */
187 
188  ///
189  /**
190  Returns a RealVect that is this RealVect with a scalar s added
191  to each component.
192  */
193  RealVect operator+ (Real s) const;
194 
195  ///
196  /**
197  Modifies this RealVect by component-wise addition by argument.
198  */
199  RealVect& operator+= (const RealVect& p);
200 
201  ///
202  /**
203  Modifies this RealVect by subtraction of a scalar from each
204  component.
205  */
207 
208  ///
209  /**
210  Modifies this RealVect by component-wise subtraction by argument.
211  */
212  RealVect& operator-= (const RealVect& p);
213 
214  ///
215  /**
216  Returns a RealVect that is this RealVect with a scalar s subtracted
217  from each component.
218  */
219  RealVect operator- (Real s) const;
220 
221  ///
222  /**
223  Modifies this RealVect by multiplying each component by a scalar.
224  */
226 
227  ///
228  /**
229 
230  */
231  Real dotProduct(const RealVect& a_rhs) const;
232 
233  ///
234  /**
235  Modifies this RealVect by component-wise multiplication by argument.
236  */
237  RealVect& operator*= (const RealVect& p);
238 
239 //XXX ///
240 //XXX /**
241 //XXX Returns component-wise product of this RealVect and argument.
242 //XXX */
243 //XXX RealVect operator* (const RealVect& p) const;
244 
245  ///
246  /**
247  Returns a RealVect that is this RealVect with each component
248  multiplied by a scalar.
249  */
250  RealVect operator* (Real s) const;
251 
252  ///
253  /**
254  Modifies this RealVect by dividing each component by a scalar.
255  */
257 
258  ///
259  /**
260  Modifies this RealVect by component-wise division by argument.
261  */
262  RealVect& operator/= (const RealVect& p);
263 
264 //XXX ///
265 //XXX /**
266 //XXX Returns component-wise quotient of this RealVect by argument.
267 //XXX */
268 //XXX RealVect operator/ (const RealVect& p) const;
269 
270  ///
271  /**
272  Returns a RealVect that is this RealVect with each component
273  divided by a scalar.
274  */
275  RealVect operator/ (Real s) const;
276 
277  ///
278  /**
279  Modifies this RealVect by multiplying each component by a scalar.
280  */
281  RealVect& scale (Real s);
282 
283  /*@}*/
284 
285  /**
286  \name Other arithmetic operators
287  */
288  /*@{*/
289 
290  ///
291  /**
292  Modifies this RealVect by taking component-wise min with RealVect
293  argument.
294  */
295  RealVect& min (const RealVect& p);
296 
297  ///
298  /**
299  Returns the RealVect that is the component-wise minimum of two
300  argument RealVects.
301  */
302  friend inline RealVect min (const RealVect& p1,
303  const RealVect& p2);
304 
305  ///
306  /**
307  Modifies this RealVect by taking component-wise max with RealVect
308  argument.
309  */
310  RealVect& max (const RealVect& p);
311 
312  ///
313  /**
314  Returns the RealVect that is the component-wise maximum of two
315  argument RealVects.
316  */
317  friend inline RealVect max (const RealVect& p1,
318  const RealVect& p2);
319 
320  /*@}*/
321 
322  /**
323  \name Unary operators
324  */
325  /*@{*/
326 
327  ///
328  /**
329  Unary plus -- for completeness.
330  */
331  RealVect operator+ () const;
332 
333  ///
334  /**
335  Unary minus -- negates all components of this RealVect.
336  */
337  RealVect operator- () const;
338 
339  ///
340  /**
341  Sum of all components of this RealVect.
342  */
343  Real sum () const;
344 
345  ///
346  /**
347  sqrt(sum squares)
348  */
349  Real vectorLength() const;
350 
351  ///
352  /**
353  sum squares--no square root
354  */
355  Real radSquared() const;
356 
357  ///
358  /**
359  Product of all components of this RealVect.
360  */
361  Real product () const;
362 
363  ///
364  /**
365  Component with the minimum value of this RealVect (returns 0 if they are all the same).
366  a_doAbs : if true then take the absolute value before comparing
367  */
368  int minDir(const bool& a_doAbs) const;
369 
370  ///
371  /**
372  Component with the maximum value of this RealVect (returns 0 if they are all the same).
373  a_doAbs : if true then take the absolute value before comparing
374  */
375  int maxDir(const bool& a_doAbs) const;
376 
377  /*@}*/
378 
379  /**
380  \name Data pointer functions
381  */
382  /*@{*/
383 
384  ///
385  /**
386  Only for sending stuff to Fortran
387  */
388  const Real* dataPtr() const;
389 
390  ///
391  /**
392  Only for sending stuff to Fortran
393  */
394  Real* dataPtr() ;
395 
396  /*@}*/
397 
398  /**
399  \name Constants
400  */
401  /*@{*/
402 
403  ///
404  /**
405  Returns a basis vector in the given coordinate direction.<br>
406  In 2-D:<br>
407  BASISREALV(0) == (1.,0.);
408  BASISREALV(1) == (0.,1.).<br>
409  In 3-D:<br>
410  BASISREALV(0) == (1.,0.,0.);
411  BASISREALV(1) == (0.,1.,0.);
412  BASISREALV(2) == (0.,0.,1.).<br>
413  Note that the coordinate directions are based at zero.
414  */
415  friend RealVect BASISREALV(int dir);
416 
417  ///
418  /**
419  This is a RealVect all of whose components are equal to zero.
420  */
421  static const RealVect Zero;
422 
423  ///
424  /**
425  This is a RealVect all of whose components are equal to one.
426  */
427  static const RealVect Unit;
428 
429  /*@}*/
430 
431  /**
432  \name Arithmetic friend functions
433  */
434  /*@{*/
435 
436  ///
437  /**
438  Returns a RealVect that is a RealVect <i>p</i> with
439  a scalar <i>s</i> added to each component.
440  */
441  friend RealVect operator+ (Real s,
442  const RealVect& p);
443 
444  ///
445  /**
446  Returns <i>s - p</i>.
447  */
448  friend RealVect operator- (Real s,
449  const RealVect& p);
450 
451  ///
452  /**
453  Returns a RealVect that is a RealVect <i>p</i> with each component
454  multiplied by a scalar <i>s</i>.
455  */
456  friend RealVect operator* (Real s,
457  const RealVect& p);
458  ///
459  /**
460  Returns a RealVect that is a RealVect <i>p</i> with each component
461  divided by a scalar <i>s</i>.
462  */
463  friend RealVect operator/ (Real s,
464  const RealVect& p);
465 
466  ///
467  /**
468  Returns component-wise sum of RealVects <i>s</i> and <i>p</i>.
469  */
470  friend RealVect operator+ (const RealVect& s,
471  const RealVect& p);
472 
473  ///
474  /**
475  Returns <i>s - p</i>.
476  */
477  friend RealVect operator- (const RealVect& s,
478  const RealVect& p);
479 
480  ///
481  /**
482  Returns component-wise product of <i>s</i> and <i>p</i>.
483  */
484  friend RealVect operator* (const RealVect& s,
485  const RealVect& p);
486  ///
487  /**
488  Returns component-wise quotient <i>p / s</i>.
489  */
490  friend RealVect operator/ (const RealVect& s,
491  const RealVect& p);
492 
493  ///
494  /**
495  Returns a RealVect obtained by multiplying each of the components
496  of the given RealVect by a scalar.
497  */
498  friend inline RealVect scale (const RealVect& p,
499  Real s);
500 
501  /*@}*/
502 
503  ///
504  /**
505  Print to the given output stream in ASCII.
506  */
507  friend std::ostream& operator<< (std::ostream& ostr,
508  const RealVect& p);
509 
510  friend class HDF5Handle;
511 
512  static size_t io_offset;
513 
514 protected:
515 
516  /**
517  The individual components of this RealVect.
518  */
520 
521 };
522 
523 #include "NamespaceFooter.H"
524 
525 #include "BaseNamespaceHeader.H"
526 
527 #include "NamespaceVar.H"
528 // template specialization of linearIn and LinearOut for RealVect
529 // RealVect spcializations of linearization
530 template < >
531 int linearSize(const CH_XDIR::RealVect& vindex);
532 
533 //VolIndex specialization of linearIn
534 template < >
535 void linearIn(CH_XDIR::RealVect& a_outputT, const void* const inBuf);
536 
537 //VolIndex specialization of linearOut
538 template < >
539 void linearOut(void* const a_outBuf, const CH_XDIR::RealVect& a_inputT);
540 
541 //Vector<RealVect> specialization
542 template < >
543 int linearSize(const Vector<CH_XDIR::RealVect>& a_input);
544 template < >
545 void linearIn(Vector<CH_XDIR::RealVect>& a_outputT, const void* const inBuf);
546 template < >
547 void linearOut(void* const a_outBuf, const Vector<CH_XDIR::RealVect>& a_inputT);
548 
549 //Vector<Vector<RealVect> > specialization
550 template < >
551 int linearSize(const Vector<Vector<CH_XDIR::RealVect> >& a_input);
552 template < >
553 void linearIn(Vector<Vector<CH_XDIR::RealVect> >& a_outputT, const void* const inBuf);
554 template < >
555 void linearOut(void* const a_outBuf, const Vector<Vector<CH_XDIR::RealVect> >& a_inputT);
556 
557 #include "BaseNamespaceFooter.H"
558 
559 #include "NamespaceHeader.H"
560 
562 {
563  CH_assert(i>=0 && i < SpaceDim);
564  return vect[i];
565 }
566 
567 inline const Real& RealVect::operator[] (int i) const
568 {
569  CH_assert(i>=0 && i < SpaceDim);
570  return vect[i];
571 }
572 
573 inline RealVect::RealVect (const RealVect &iv)
574 {
575  D_EXPR6(vect[0]=iv.vect[0], vect[1]=iv.vect[1], vect[2]=iv.vect[2],
576  vect[3]=iv.vect[3], vect[4]=iv.vect[4], vect[5]=iv.vect[5]);
577 }
578 
579 inline
580 RealVect&
582 {
583  D_EXPR6(vect[0] -= s, vect[1] -= s, vect[2] -= s,
584  vect[3] -= s, vect[4] -= s, vect[5] -= s);
585  return *this;
586 }
587 
588 inline
589 RealVect&
591 {
592  D_EXPR6(vect[0] -= p[0], vect[1] -= p[1], vect[2] -= p[2],
593  vect[3] -= p[3], vect[4] -= p[4], vect[5] -= p[5]);
594  return *this;
595 }
596 
597 inline
598 RealVect
600 {
601  return RealVect(*this);
602 }
603 
604 inline
605 RealVect
607 {
608  return RealVect(D_DECL6(-vect[0], -vect[1], -vect[2],
609  -vect[3], -vect[4], -vect[5] ));
610 }
611 
612 inline
613 RealVect&
615 {
616  D_EXPR6(vect[0] *= s, vect[1] *= s, vect[2] *= s,
617  vect[3] *= s, vect[4] *= s, vect[5] *= s);
618  return *this;
619 }
620 
621 inline
622 Real
624 {
625  return D_TERM6(vect[0], + vect[1], + vect[2], +
626  vect[3], + vect[4], + vect[5]);
627 }
628 
629 inline
630 Real
632 {
633  Real len = this->radSquared();
634  len = sqrt(len);
635 
636  return len;
637 }
638 
639 inline
640 Real
642 {
643  Real len = 0;
644  for (int idir = 0; idir < SpaceDim; idir++)
645  {
646  len = len + vect[idir]*vect[idir];
647  }
648 
649  return len;
650 }
651 
652 inline
653 Real
655 {
656  return D_TERM6(vect[0], * vect[1], * vect[2], *
657  vect[3], * vect[4], * vect[5]);
658 }
659 
660 inline
661 RealVect
662 scale (const RealVect& p,
663  Real s)
664 {
665  return RealVect(D_DECL6(s * p[0], s * p[1], s * p[2],
666  s * p[3], s * p[4], s * p[5]));
667 }
668 
669 inline
670 bool
672 {
673  return D_TERM6(vect[0] < p[0], && vect[1] < p[1], && vect[2] < p[2],
674  && vect[3] < p[3], && vect[4] < p[4], && vect[5] < p[5]);
675 }
676 
677 inline
678 bool
680 {
681  return D_TERM6(vect[0] <= p[0], && vect[1] <= p[1], && vect[2] <= p[2],
682  && vect[3] <= p[3], && vect[4] <= p[4], && vect[5] <= p[5]);
683 }
684 
685 inline
686 bool
688 {
689  return D_TERM6(vect[0] > p[0], && vect[1] > p[1], && vect[2] > p[2],
690  && vect[3] > p[3], && vect[4] > p[4], && vect[5] > p[5]);
691 }
692 
693 inline
694 bool
696 {
697  return D_TERM6(vect[0] >= p[0], && vect[1] >= p[1], && vect[2] >= p[2],
698  && vect[3] >= p[3], && vect[4] >= p[4], && vect[5] >= p[5]);
699 }
700 
701 inline
702 RealVect&
704 {
705  D_EXPR6(vect[0] = Min(vect[0], p.vect[0]),
706  vect[1] = Min(vect[1], p.vect[1]),
707  vect[2] = Min(vect[2], p.vect[2]),
708  vect[3] = Min(vect[3], p.vect[3]),
709  vect[4] = Min(vect[4], p.vect[4]),
710  vect[5] = Min(vect[5], p.vect[5]));
711  return *this;
712 }
713 
714 inline
715 RealVect&
717 {
718  D_EXPR6(vect[0] = Max(vect[0], p.vect[0]),
719  vect[1] = Max(vect[1], p.vect[1]),
720  vect[2] = Max(vect[2], p.vect[2]),
721  vect[3] = Max(vect[3], p.vect[3]),
722  vect[4] = Max(vect[4], p.vect[4]),
723  vect[5] = Max(vect[5], p.vect[5]));
724  return *this;
725 }
726 
727 inline
728 RealVect
729 min (const RealVect& p1,
730  const RealVect& p2)
731 {
732  RealVect p(p1);
733  return p.min(p2);
734 }
735 
736 inline
737 RealVect
738 max (const RealVect& p1,
739  const RealVect& p2)
740 {
741  RealVect p(p1);
742  return p.max(p2);
743 }
744 
745 extern RealVect BASISREALV(int idir);
746 
747 #include "NamespaceFooter.H"
748 
749 #endif
#define D_DECL6(a, b, c, d, e, f)
Definition: CHArray.H:39
RealVect & operator/=(Real s)
RealVect operator-() const
Definition: RealVect.H:606
RealVect & operator=(const RealVect &rhs)
RealVect operator*(Real s) const
#define D_TERM6(a, b, c, d, e, f)
Definition: CHArray.H:40
#define CH_assert(cond)
Definition: CHArray.H:37
Real radSquared() const
Definition: RealVect.H:641
RealVect(const IntVect &iv)
Definition: RealVect.H:84
int maxDir(const bool &a_doAbs) const
RealVect & scale(Real s)
Definition: RealVect.H:614
const Real * dataPtr() const
bool operator>=(const RealVect &p) const
Definition: RealVect.H:695
RealVect operator+() const
Definition: RealVect.H:599
Real sum() const
Definition: RealVect.H:623
Real & operator[](int i)
Definition: RealVect.H:561
void linearIn(CH_XDIR::RealVect &a_outputT, const void *const inBuf)
bool operator!=(const RealVect &p) const
Real product() const
Definition: RealVect.H:654
const int SpaceDim
Definition: SPACE.H:38
Definition: IndexTM.H:36
static const RealVect Unit
Definition: RealVect.H:427
Real dotProduct(const RealVect &a_rhs) const
Real vectorLength() const
Definition: RealVect.H:631
static const RealVect Zero
Definition: RealVect.H:421
bool operator<(const RealVect &p) const
Definition: RealVect.H:671
RealVect & max(const RealVect &p)
Definition: RealVect.H:716
double Real
Definition: REAL.H:33
void linearOut(void *const a_outBuf, const CH_XDIR::RealVect &a_inputT)
RealVect & operator*=(Real s)
bool operator<=(const RealVect &p) const
Definition: RealVect.H:679
bool operator>(const RealVect &p) const
Definition: RealVect.H:687
friend std::ostream & operator<<(std::ostream &ostr, const RealVect &p)
A Real vector in SpaceDim-dimensional space.
Definition: RealVect.H:41
RealVect & min(const RealVect &p)
Definition: RealVect.H:703
Handle to a particular group in an HDF file.
Definition: CH_HDF5.H:294
T Min(const T &a_a, const T &a_b)
Definition: Misc.H:26
An integer Vector in SpaceDim-dimensional space.
Definition: CHArray.H:42
RealVect operator/(Real s) const
int minDir(const bool &a_doAbs) const
T Max(const T &a_a, const T &a_b)
Definition: Misc.H:39
RealVect & operator+=(Real s)
Real vect[SpaceDim]
Definition: RealVect.H:519
bool operator==(const RealVect &p) const
int linearSize(const CH_XDIR::RealVect &vindex)
RealVect & operator-=(Real s)
Definition: RealVect.H:581
friend RealVect BASISREALV(int dir)
static size_t io_offset
Definition: RealVect.H:512