Carmen Code
 All Classes Files Functions Variables Macros Pages
Functions
TimeEvolution.cpp File Reference

Time evolution for finite volume with multiresolution. More...

#include "Carmen.h"
Include dependency graph for TimeEvolution.cpp:

Functions

void TimeEvolution (Node *Root)
 Computes a time evolution on the tree structure, the root node being Root. Only for multiresolution computations. More...
 
void TimeEvolution (FineMesh *Root)
 Computes a time evolution on the regular fine mesh Root. Only for finite volume computations. More...
 

Detailed Description

Time evolution for finite volume with multiresolution.

Function Documentation

void TimeEvolution ( Node Root)

Computes a time evolution on the tree structure, the root node being Root. Only for multiresolution computations.

Parameters
RootRoot node
Returns
void
20 {
21 
22  // --- Smooth data ---
23 
24  if ( SmoothCoeff != 0.)
25  Root->smooth();
26 
27  // --- Store cell-average values of leaves ---
28  Root->store();
29 
30  // --- Refresh tree structure ---
31  RefreshTree(Root);
32 
33  for (StepNo = 1; StepNo <= StepNb; StepNo++)
34  {
35  // --- Compute divergence ---
36  Root->computeDivergence();
37  // --- Runge-Kutta step ---
38  Root->RungeKutta();
39  // --- Divergence cleaning source-terms
40  Root->computeCorrection();
41 
42  }
43 
44  // --- Refresh tree structure ---
45  RefreshTree(Root);
46 
47 
48 
49  // --- Check stability ---
50  Root->checkStability();
51 
52  // --- Compute integral values ---
53  Root->computeIntegral();
54 
55  // --- Compute total number of cells and leaves ---
56 
59  //cout<<"eigen= "<<Eigenvalue<<endl;
60  // --- Compute elapsed time and adapt time step ---
64 
65  // --- Compute divergence-free correction constant
66  //ch = CFL*SpaceStep/TimeStep;
68 
69 }
void RefreshTree(Node *Root)
Refresh the tree structure, i.e. compute the cell-averages of the internal nodes by projection and th...
Definition: RefreshTree.cpp:22
real ch
Definition: Parameters.cpp:140
real EigenvalueZ
Definition: Parameters.cpp:160
int StepNo
Definition: Parameters.cpp:167
void RungeKutta()
Computes one Runge-Kutta step.
Definition: Node.cpp:2345
int LeafNb
Definition: Parameters.cpp:175
real Eigenvalue
Definition: Parameters.cpp:157
real TotalLeafNb
Definition: Parameters.cpp:176
real SmoothCoeff
Definition: Parameters.cpp:99
real SpaceStep
Definition: Parameters.cpp:156
void computeCorrection()
Computes velocity gradient (only for Navier-Stokes).
Definition: Node.cpp:2188
real ElapsedTime
Definition: Parameters.cpp:95
int StepNb
Definition: Parameters.cpp:36
int CellNb
Definition: Parameters.cpp:174
real TotalCellNb
Definition: Parameters.cpp:177
void checkStability()
Checks if the computation is numerically unstable, i.e. if one of the cell-averages is overflow...
Definition: Node.cpp:2430
real TimeStep
Definition: Parameters.cpp:40
void AdaptTimeStep()
Adapts time step when required.
Definition: AdaptTimeStep.cpp:24
void computeIntegral()
Computes integral values like e.g. flame velocity, global error, etc.
Definition: Node.cpp:2470
bool ConstantTimeStep
Definition: Parameters.cpp:94
real EigenvalueY
Definition: Parameters.cpp:159
void computeDivergence()
Computes the divergence vector with the space discretization scheme.
Definition: Node.cpp:1977
void smooth()
Deletes the details in the highest level.
Definition: Node.cpp:2874
real EigenvalueX
Definition: Parameters.cpp:158
#define Max(x, y)
Definition: Carmen.h:54
real CFL
Definition: Parameters.cpp:39
void store()
Stores cell-average values into temporary cell-average values.
Definition: Node.cpp:1913
void TimeEvolution ( FineMesh Root)

Computes a time evolution on the regular fine mesh Root. Only for finite volume computations.

Parameters
RootFine mesh
Returns
void
78 {
79 
80  // --- Store cell-average values into temporary ---
81  Root->store();
82 
83  for (StepNo = 1; StepNo <= StepNb; StepNo++)
84  {
85  // --- Compute divergence for neighbour cells ---
86  //The same conception with background computations, see upper...
87  Root->computeDivergence(1);
88  // --- Runge-Kutta step for neighbour cells ---
89  Root->RungeKutta(1);
90  // --- Divergence cleaning source term
91  Root->computeCorrection(1);
92  // --- Start inter-CPU exchanges ---
93  CPUExchange(Root, SendQ);
94  // --- Compute divergence for internal cells ---
95  Root->computeDivergence(0);
96  // --- Runge-Kutta step for internal cells ---
97  Root->RungeKutta(0);
98  // --- Divergence cleaning source term
99  Root->computeCorrection(0);
100 
101 #if defined PARMPI
102  CommTimer.start(); //Communication Timer Start
103  //Waiting while inter-CPU exchanges are finished
104  if (MPIRecvType == 1) //for nonblocking recive...
105  MPI_Waitall(4*Dimension,req,st);
106  CommTimer.stop();
107 #endif
108 
109  }
110 
111  // --- Check stability ---
112  Root->checkStability();
113 
114  // --- Compute integral values ---
115  Root->computeIntegral();
116 
117  // --- Compute elapsed time and adapt time step ---
118 
119  if (!ComputeCPUTimeRef)
120  {
124  // --- Compute divergence-free correction constant
125  //ch = CFL*SpaceStep/TimeStep;
127  }
128 
129  // --- Compute time-average values ---
130 
131  if (TimeAveraging)
132  Root->computeTimeAverage();
133 
134 }
void computeCorrection(int)
Computes divergence cleaning source term (only for MHD).
Definition: FineMesh.cpp:540
void checkStability() const
Checks if the computation is numerically unstable, i.e. if one of the cell-averages is overflow...
Definition: FineMesh.cpp:927
real ch
Definition: Parameters.cpp:140
real EigenvalueZ
Definition: Parameters.cpp:160
bool TimeAveraging
Definition: Parameters.cpp:60
bool ComputeCPUTimeRef
Definition: Parameters.cpp:62
int MPIRecvType
Definition: Parameters.cpp:239
int StepNo
Definition: Parameters.cpp:167
double stop()
Stop timer and, if asked, returns CPU time from previous start in seconds.
Definition: Timer.cpp:77
void computeTimeAverage()
Computes the time-average value in every cell.
Definition: FineMesh.cpp:1803
real Eigenvalue
Definition: Parameters.cpp:157
void computeIntegral()
Computes integral values like e.g. flame velocity, global error, etc.
Definition: FineMesh.cpp:965
real SpaceStep
Definition: Parameters.cpp:156
int Dimension
Definition: Parameters.cpp:74
real ElapsedTime
Definition: Parameters.cpp:95
int StepNb
Definition: Parameters.cpp:36
Timer CommTimer
Definition: Parameters.cpp:241
void RungeKutta(int)
Computes the Runge-Kutta step.
Definition: FineMesh.cpp:807
void start()
Starts timer.
Definition: Timer.cpp:62
real TimeStep
Definition: Parameters.cpp:40
void AdaptTimeStep()
Adapts time step when required.
Definition: AdaptTimeStep.cpp:24
bool ConstantTimeStep
Definition: Parameters.cpp:94
real EigenvalueY
Definition: Parameters.cpp:159
void store()
Stores cell-average values into temporary cell-average values.
Definition: FineMesh.cpp:882
void computeDivergence(int)
Computes the divergence vector with the space discretization scheme.
Definition: FineMesh.cpp:426
int SendQ
Definition: Parameters.cpp:252
real EigenvalueX
Definition: Parameters.cpp:158
#define Max(x, y)
Definition: Carmen.h:54
real CFL
Definition: Parameters.cpp:39
void CPUExchange(FineMesh *Root, int)
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:350

Here is the caller graph for this function: