Carmen Code
 All Classes Files Functions Variables Macros Pages
Node.h
Go to the documentation of this file.
1 /***************************************************************************
2  Node.h - description
3  -------------------
4  begin : Thu Jun 7 2001
5  copyright : (C) 2001 by Olivier Roussel & Alexei Tsigulin
6  email : roussel@ict.uni-karlsruhe.de, lpsoft@mail.ru
7  ***************************************************************************/
8 
9 /***************************************************************************
10  * *
11  * This program is free software; you can redistribute it and/or modify *
12  * it under the terms of the GNU General Public License as published by *
13  * the Free Software Foundation; either version 2 of the License, or *
14  * (at your option) any later version. *
15  * *
16  ***************************************************************************/
17 
18 #include "PreProcessor.h"
19 
38 class Node
39 {
40 /*
41 ______________________________________________________________________________________________
42 
43 PUBLIC METHODS
44 
45  Constructor and distructor
46 ______________________________________________________________________________________________
47 
48 */
49 
50 public:
51 
64  Node(const int l=0, const int i=0, const int j=0, const int k=0);
65 
72  ~Node();
73 /*
74 ______________________________________________________________________________________________
75 
76  Get procedures
77 ______________________________________________________________________________________________
78 
79 */
80 
86  inline int cells() const;
87 
93  inline int leaves() const;
94 
95 /*
96 ______________________________________________________________________________________________
97 
98  Multiresolution procedures
99 ______________________________________________________________________________________________
100 
101 */
102 
109  int adapt();
110 
117  void checkGradedTree();
118 
119 
125  void initValue();
126 
132  void addLevel();
133 
141  Cell* project();
142 
150  void fillVirtualChildren();
151 /*
152 ______________________________________________________________________________________________
153 
154  Time evolution procedures
155 ______________________________________________________________________________________________
156 
157 */
158 
164 void store();
165 
171 void storeGrad();
172 
178 void computeDivergence();
179 
185 void RungeKutta();
186 
192 void computeIntegral();
193 
199 void computeGradient();
200 
206 void computeCorrection();
207 
215 void checkStability();
216 /*
217 ______________________________________________________________________________________________
218 
219  Output procedures
220 ______________________________________________________________________________________________
221 
222 */
223 
230 void writeTree(const char* FileName) const;
231 
239 void writeAverage(const char* FileName);
240 
247 void writeMesh(const char* FileName) const;
248 
255 void writeHeader(const char* FileName) const;
256 
264 void writeFineGrid(const char* FileName, const int L=ScaleNb) const;
265 
266 /*
267 ______________________________________________________________________________________________
268 
269  Backup-restore procedure (to restart a computation)
270 ______________________________________________________________________________________________
271 
272 */
279 void backup();
280 
287 void restore();
288 
295 void restoreFineMesh();
296 
302 void smooth();
303 
304 /*
305 ______________________________________________________________________________________________
306 
307 PRIVATE METHODS
308 ______________________________________________________________________________________________
309 
310  Multiresolution procedures
311 ______________________________________________________________________________________________
312 
313 */
314 
315 private:
324  void split(const bool init=false);
325 
331 void splitAll();
332 
338 void combine();
339 
345 Vector predict() const;
346 
353 Vector predictTempAverage() const;
354 
360 Matrix predictGradient() const;
361 
367 inline void computeTimeInterpolation();
368 
374 inline void computeTimeExtrapolation();
375 
381 inline void storeTimeEvolution();
382 
388 inline void setInternalNode();
389 
395 inline void setSimpleLeaf();
396 
402 inline void setLeafWithVirtualChildren();
403 
409 inline void setVirtualLeaf();
410 
416 inline bool isInternalNode() const;
417 
423 inline bool isLeaf() const;
424 
430 inline bool isSimpleLeaf() const;
431 
437 inline bool isLeafWithVirtualChildren() const;
438 
444 inline bool isVirtualLeaf() const;
445 
451 inline bool hasChildren() const;
452 
458 inline bool isParentOfLeaf() const;
459 
465 bool isInsideBoundary() const;
466 
472 bool isInFluid() const;
473 
479 bool isOnBoundary() const;
480 
487 inline bool isBeginTimeCycle() const;
488 
495 inline bool isEndTimeCycle() const;
496 
502 inline bool requiresTimeEvolution() const;
503 
509 inline bool requiresTimeInterpolation() const;
510 
516 inline bool requiresDivergenceComputation() const;
517 
523 inline bool detailIsSmall() const;
524 /*
525 ______________________________________________________________________________________________
526 
527  Get procedures
528 ______________________________________________________________________________________________
529 
530 */
531 
543 Node* node(int l, int i, int j = 0, int k = 0) const;
544 
554 Cell* cell(int l, int i, int j = 0, int k = 0) const;
555 
561 inline Node* parent() const;
562 
568 inline Cell* parentCell() const;
569 
578 inline Node* uncle(const int i, const int j=0, const int k=0) const;
579 
588 inline Cell* uncleCell(const int i, const int j=0, const int k=0) const;
589 
598 inline Node* cousin(const int i, const int j=0, const int k=0) const;
599 
608 inline Cell* cousinCell(const int i, const int j=0, const int k=0) const;
609 
618 inline Node* child(const int i, const int j=0, const int k=0) const;
619 
628 inline Cell* childCell(const int i, const int j=0, const int k=0) const;
629 
630 /*
631 ______________________________________________________________________________________________
632 
633  Procedures on virtual children
634 ______________________________________________________________________________________________
635 
636 */
637 
643 void deleteVirtualChildren();
644 
653 void makeVirtualChildren(bool init=false);
654 
655 /*
656 ______________________________________________________________________________________________
657 
658 PRIVATE VARIABLES
659 ______________________________________________________________________________________________
660 
661 */
662  static Node *Root;
663 // static Node Boundary; /*!< Boundary cell*/
664 
665  int Nl;
666  int Ni, Nj, Nk;
668  Node **Child;
669  Cell ThisCell;
671  byte Flag;
672 };
673 
674 /*
675 ______________________________________________________________________________________________
676 
677  INLINE FUNCTION
678 ______________________________________________________________________________________________
679 
680  Returns the number of cells
681 ______________________________________________________________________________________________
682 
683 */
684 inline int Node::cells() const
685 {
686  return CellNb;
687 }
688 /*
689 ______________________________________________________________________________________________
690 
691  Get number of leaves
692 ______________________________________________________________________________________________
693 
694 */
695 inline int Node::leaves() const
696 {
697  return LeafNb;
698 }
699 
700 /*
701 ______________________________________________________________________________________________
702 
703  Returns pointer to child cell using relative position
704 ______________________________________________________________________________________________
705 
706 */
707 
708 inline Cell* Node::childCell(const int i, const int j, const int k) const
709 {
710  return cell(Nl+1, 2*Ni+i, 2*Nj+j, 2*Nk+k);
711 }
712 
713 /*
714 ______________________________________________________________________________________________
715 
716  Returns pointer to child node using relative position
717 ______________________________________________________________________________________________
718 
719 */
720 
721 inline Node* Node::child(const int i, const int j, const int k) const
722 {
723  return node(Nl+1, 2*Ni+i, 2*Nj+j, 2*Nk+k);
724 }
725 /*
726 ______________________________________________________________________________________________
727 
728  This function computes the interpolation in time in the node
729 ______________________________________________________________________________________________
730 
731 */
732 
733 inline void Node::computeTimeInterpolation()
734 {
735  ThisCell.setAverage(ThisCell.tempAverage() + ThisCell.average());
736 }
737 
738 /*
739 ______________________________________________________________________________________________
740 
741  Returns pointer to uncle cell using relative position
742 ______________________________________________________________________________________________
743 
744 */
745 
746 inline Cell* Node::cousinCell(const int i, const int j, const int k) const
747 {
748  return cell(Nl, Ni+i, Nj+j, Nk+k);
749 }
750 
751 /*
752 ______________________________________________________________________________________________
753 
754  Returns pointer to uncle cell using relative position
755 ______________________________________________________________________________________________
756 
757 */
758 
759 inline Node* Node::cousin(const int i, const int j, const int k) const
760 {
761  return node(Nl, Ni+i, Nj+j, Nk+k);
762 }
763 
764 /*
765 ______________________________________________________________________________________________
766 
767  Return true if node is has children (real or virtual)
768 ______________________________________________________________________________________________
769 
770 */
771 inline bool Node::hasChildren() const
772 {
773  return (Flag==0 || Flag==2);
774 }
775 
776 /*
777 ______________________________________________________________________________________________
778 
779  Sets the current node to internal node
780 ______________________________________________________________________________________________
781 
782 */
783 inline void Node::setInternalNode()
784 {
785  Flag = 0;
786  return;
787 }
788 
789 /*
790 ______________________________________________________________________________________________
791 
792  Sets the current node to simple leaf
793 ______________________________________________________________________________________________
794 
795 */
796 inline void Node::setSimpleLeaf()
797 {
798  Flag = 1;
799  return;
800 }
801 /*
802 ______________________________________________________________________________________________
803 
804  Sets the current node to leaf with virtual children
805 ______________________________________________________________________________________________
806 
807 */
808 inline void Node::setLeafWithVirtualChildren()
809 {
810  Flag = 2;
811  return;
812 }
813 /*
814 ______________________________________________________________________________________________
815 
816  Sets the current node to virtual leaf
817 ______________________________________________________________________________________________
818 
819 */
820 inline void Node::setVirtualLeaf()
821 {
822  Flag = 3;
823  return;
824 }
825 /*
826 ______________________________________________________________________________________________
827 
828 Returns true if the leaf is at the begining of a time cycle.
829 Only useful when TimeAdaptivity = true
830 ______________________________________________________________________________________________
831 
832 */
833 inline bool Node::isBeginTimeCycle() const
834 {
835  return ((IterationNo-1)%(1<<(TimeAdaptivityFactor*(ScaleNb-Nl))) == 0);
836 }
837 
838 /*
839 ______________________________________________________________________________________________
840 
841 Returns true if the leaf is at the end of a time cycle.
842 Only useful when TimeAdaptivity = true
843 ______________________________________________________________________________________________
844 
845 */
846 inline bool Node::isEndTimeCycle() const
847 {
848  return (IterationNo%(1<<(TimeAdaptivityFactor*(ScaleNb-Nl))) == 0);
849 }
850 
851 /*
852 ______________________________________________________________________________________________
853 
854  Returns pointer to parent node
855 ______________________________________________________________________________________________
856 
857 */
858 inline Node* Node::parent() const
859 {
860  return node(Nl-1, (Ni+4)/2-2, (Nj+4)/2-2, (Nk+4)/2-2);
861 }
862 
863 /*
864 ______________________________________________________________________________________________
865 
866  Returns pointer to parent cell
867 ______________________________________________________________________________________________
868 
869 */
870 inline Cell* Node::parentCell() const
871 {
872  return cell(Nl-1, (Ni+4)/2-2, (Nj+4)/2-2, (Nk+4)/2-2);
873 }
874 
875 /*
876 ______________________________________________________________________________________________
877 
878  Returns pointer to uncle node using relative position
879 ______________________________________________________________________________________________
880 
881 */
882 inline Node* Node::uncle(const int i, const int j, const int k) const
883 {
884  return node(Nl-1, (Ni+4)/2-2 + i, (Nj+4)/2-2 + j, (Nk+4)/2-2 + k);
885 }
886 
887 /*
888 ______________________________________________________________________________________________
889 
890  Returns pointer to uncle cell using relative position
891 ______________________________________________________________________________________________
892 
893 */
894 inline Cell* Node::uncleCell(const int i, const int j, const int k) const
895 {
896  return cell(Nl-1, (Ni+4)/2-2 + i, (Nj+4)/2-2 + j, (Nk+4)/2-2 + k);
897 }
898 
899 /*
900 ______________________________________________________________________________________________
901 
902  Returns true if node is an internal node (not a leaf)
903 ______________________________________________________________________________________________
904 
905 */
906 inline bool Node::isInternalNode() const
907 {
908  return (Flag==0);
909 }
910 
911 /*
912 ______________________________________________________________________________________________
913 
914  Returns true if node is a leaf (with or without virtual children)
915 ______________________________________________________________________________________________
916 
917 */
918 inline bool Node::isLeaf() const
919 {
920  return (Flag==1 || Flag==2);
921 }
922 
923 /*
924 ______________________________________________________________________________________________
925 
926  Returns true if node is a simple leaf (without virtual children)
927 ______________________________________________________________________________________________
928 
929 */
930 inline bool Node::isSimpleLeaf() const
931 {
932  return (Flag==1);
933 }
934 
935 /*
936 ______________________________________________________________________________________________
937 
938  Returns true if node is a leaf with virtual children
939 ______________________________________________________________________________________________
940 
941 */
942 inline bool Node::isLeafWithVirtualChildren() const
943 {
944  return (Flag==2);
945 }
946 
947 /*
948 ______________________________________________________________________________________________
949 
950  Return true if node is a virtual leaf
951 ______________________________________________________________________________________________
952 
953 */
954 inline bool Node::isVirtualLeaf() const
955 {
956  return (Flag==3);
957 }
958 
959 /*
960 ______________________________________________________________________________________________
961 
962  This function stores the result of the time evolution in Qlow, and interpolates
963  the internediary state between [n+1] and [n], which is stored in Q.
964 ______________________________________________________________________________________________
965 
966 */
967 
968 inline void Node::storeTimeEvolution()
969 {
970 
971  // Qlow <- Q
972  ThisCell.setLowAverage( ThisCell.average());
973 
974  // Q <- 1/2 (Q + Qs)
975  ThisCell.setAverage( 0.5*(ThisCell.average() + ThisCell.tempAverage()) );
976 }
977 
978 /*
979 ______________________________________________________________________________________________
980 
981  Returns the extrapolation in time. Only for local time stepping
982 ______________________________________________________________________________________________
983 
984 */
985 
986 inline void Node::computeTimeExtrapolation()
987 {
988  Vector A(QuantityNb);
989 
990  A = ThisCell.average();
991  ThisCell.setAverage(ThisCell.average() + ThisCell.average() - ThisCell.tempAverage());
992  ThisCell.setTempAverage(A);
993 }
994 
995 /*
996 ______________________________________________________________________________________________
997 
998  This function return "true" if this node requires an interpolation in time, instead of a full time evolution
999 ______________________________________________________________________________________________
1000 
1001 */
1002 
1003 inline bool Node::requiresDivergenceComputation() const
1004 {
1005  if (UseBoundaryRegions)
1006  if (isInsideBoundary()) return false;
1007 
1008  if (TimeAdaptivity)
1009  // With time adaptivity, compute divergence at the begining and the end of time cycles
1010  return (isLeaf() && (isBeginTimeCycle()|| isEndTimeCycle()) );
1011  else
1012  // Without time adaptivity, perform TimeEvolution on leaves whatever IterationNo
1013  return (isLeaf()) ;
1014 }
1015 
1016 /*
1017 ______________________________________________________________________________________________
1018 
1019  This function return "true" if this node requires a time evolution procedure (no interpolation)
1020 ______________________________________________________________________________________________
1021 
1022 */
1023 inline bool Node::requiresTimeEvolution() const
1024 {
1025  if (UseBoundaryRegions)
1026  if (isInsideBoundary()) return false;
1027 
1028  if (TimeAdaptivity)
1029  // With time adaptivity, perform TimeEvolution on leaves every 2^a(L-l) iterations
1030  return (isLeaf() && isBeginTimeCycle());
1031  else
1032  // Without time adaptivity, perform TimeEvolution on leaves whatever IterationNo
1033  return (isLeaf()) ;
1034 }
1035 
1036 /*
1037 ______________________________________________________________________________________________
1038 
1039  This function return "true" if this node requires an interpolation in time, instead of a full time evolution
1040 ______________________________________________________________________________________________
1041 
1042 */
1043 inline bool Node::requiresTimeInterpolation() const
1044 {
1045  if (UseBoundaryRegions)
1046  if (isInsideBoundary()) return false;
1047 
1048  // When the time adaptivity is not used, never perform interpolation in time
1049  if (!TimeAdaptivity) return false;
1050 
1051  // With time adaptivity, permform interpolation on leaves where no time evolution is made
1052  return (isLeaf() && isEndTimeCycle());
1053 }
1054 
1055 
1056 
1057 
1058 
void storeGrad()
Stores gradient values into temporary gradient values.
Definition: Node.cpp:1950
void writeTree(const char *FileName) const
Writes tree structure into file FileName. Only for debugging.
Definition: Node.cpp:955
bool TimeAdaptivity
Definition: Parameters.cpp:85
An object Cell contains all the informations of a cell for both multiresolution and finite volume com...
Definition: Cell.h:41
int ScaleNb
Definition: Parameters.cpp:87
An object Node is an element of a graded tree structure, used for multiresolution computations...
Definition: Node.h:38
int QuantityNb
Definition: Parameters.cpp:171
void computeGradient()
Computes velocity gradient (only for Navier-Stokes).
Definition: Node.cpp:2266
Standard class for every vector in Carmen.
Definition: Vector.h:29
void writeFineGrid(const char *FileName, const int L=ScaleNb) const
Writes cell-average values on a regular grid of level L into file FileName.
Definition: Node.cpp:1174
void RungeKutta()
Computes one Runge-Kutta step.
Definition: Node.cpp:2345
int LeafNb
Definition: Parameters.cpp:175
int leaves() const
Returns the number of leaves in the tree.
Definition: Node.h:695
void backup()
Backs up the tree structure and the cell-averages into a file carmen.bak. In further computations...
Definition: Node.cpp:2716
void restoreFineMesh()
Restores the tree structure and the cell-averages from the file carmen.bak in FineMesh format...
Definition: Node.cpp:2821
void checkGradedTree()
Checks if the tree is graded. If not, an error is emitted. Only for debugging.
Definition: Node.cpp:2630
int IterationNo
Definition: Parameters.cpp:168
void computeCorrection()
Computes velocity gradient (only for Navier-Stokes).
Definition: Node.cpp:2188
#define byte
Definition: PreProcessor.h:45
void writeMesh(const char *FileName) const
Writes mesh data for Gnuplot into file FileName.
Definition: Node.cpp:1819
void initValue()
Computes the initial value.
Definition: Node.cpp:124
void fillVirtualChildren()
Fills the cell-average values of every virtual leaf with values predicted from its parent and uncles...
Definition: Node.cpp:624
int adapt()
Computes the details in the leaves and its parent nodes and, in function of the threshold Tolerance...
Definition: Node.cpp:691
int CellNb
Definition: Parameters.cpp:174
void setTempAverage(const int QuantityNo, const real UserAverage)
Identical to setAverage (int QuantityNo, real UserAverage), but for the vector of the temporary cell-...
Definition: Cell.h:945
void writeAverage(const char *FileName)
Writes cell-average values in multiresolution representation and the corresponding mesh into file Fil...
Definition: Node.cpp:1089
real average(const int QuantityNo) const
Returns the component no. QuantityNo of the cell-average value.
Definition: Cell.h:1128
Standard class for every matrix in Carmen.
Definition: Matrix.h:28
void setLowAverage(const int QuantityNo, const real UserAverage)
Identical to setAverage (int QuantityNo, real UserAverage), but for the vector of the cell-average va...
Definition: Cell.h:969
real tempAverage(const int QuantityNo) const
Returns the component no. QuantityNo of the temporary cell-average value.
Definition: Cell.h:1144
void checkStability()
Checks if the computation is numerically unstable, i.e. if one of the cell-averages is overflow...
Definition: Node.cpp:2430
void computeIntegral()
Computes integral values like e.g. flame velocity, global error, etc.
Definition: Node.cpp:2470
int cells() const
Returns the number of cells in the tree.
Definition: Node.h:684
~Node()
Distructor of Node class. Removes the node from the tree structure. If the node is not a leaf...
Definition: Node.cpp:95
void computeDivergence()
Computes the divergence vector with the space discretization scheme.
Definition: Node.cpp:1977
bool UseBoundaryRegions
Definition: Parameters.cpp:80
int TimeAdaptivityFactor
Definition: Parameters.cpp:86
void writeHeader(const char *FileName) const
Writes header for Data Explorer into file FileName.
Definition: Node.cpp:1028
void addLevel()
Adds levels when needed.
Definition: Node.cpp:191
void smooth()
Deletes the details in the highest level.
Definition: Node.cpp:2874
Node(const int l=0, const int i=0, const int j=0, const int k=0)
Constructor of Node class. Generates a new node at the position (l, i, j, k) in the tree structure...
Definition: Node.cpp:38
void restore()
Restores the tree structure and the cell-averages from the file carmen.bak. This file was created by ...
Definition: Node.cpp:2761
void setAverage(const int QuantityNo, const real UserAverage)
Sets the cell-average of the quantity QuantityNo to UserAverage. Example:
Definition: Cell.h:921
Cell * project()
Computes the cell-average values of all nodes that are not leaves by projection from the cell-average...
Definition: Node.cpp:661
void store()
Stores cell-average values into temporary cell-average values.
Definition: Node.cpp:1913