Carmen Code
 All Classes Files Functions Variables Macros Pages
Carmen.h
Go to the documentation of this file.
1 /***************************************************************************
2  Carmen.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 
22 #include "carmen.pre"
23 
24 
25 #if defined PARMPI
26 #include <mpi.h>
27 #endif
28 
29 #include "PreProcessor.h"
30 
31 
32 #include <stdio.h>
33 #include <stdlib.h>
34 #include <iostream>
35 #include <math.h>
36 
37 using namespace std;
38 
39 
40 #include "Vector.h"
41 #include "Matrix.h"
42 #include "Timer.h"
43 #include "Parameters.h"
44 #include "PrintGrid.h"
45 #include "TimeAverageGrid.h"
46 #include "Cell.h"
47 #include "Node.h"
48 #include "FineMesh.h"
49 
50 
54 #define Max(x,y) (((x) > (y)) ? (x):(y))
55 
58 #define Max3(x,y,z) (Max((x),Max((y),(z))))
59 
62 #define Min(x,y) (((x) < (y)) ? (x):(y))
63 
66 #define Min3(x,y,z) (Min((x),Min((y),(z))))
67 
70 #define power2(x) ((x)*(x))
71 
74 #define power3(x) ((x)*(x)*(x))
75 
78 #define Abs(x) ( ((x) < 0)? -(x):(x) )
79 
80 /*
81 ______________________________________________________________________________________________
82 
83 External functions
84 ______________________________________________________________________________________________
85 
86 */
87 
93 void AdaptTimeStep();
94 
107 int BC(int i, int AxisNo, int N=(1<<ScaleNb));
108 
122 int BoundaryRegion(const Vector& X);
123 
131 void Backup(Node* Root);
132 
140 void Backup(FineMesh* Root);
141 
150 double CPUTimeRef(int iterations, int scales);
151 
159 real ComputedTolerance(const int ScaleNo);
160 
167 int DigitNumber(int arg);
168 
179 int FileWrite(FILE *f, const char* format, real arg);
180 
194 Vector Flux(Cell& Cell1, Cell& Cell2, Cell& Cell3, Cell& Cell4, int AxisNo);
195 
211 void GetBoundaryCells(Cell& Cell1, Cell& Cell2, Cell& Cell3, Cell& Cell4,
212  Cell& C1, Cell& C2, Cell& C3, Cell& C4, int AxisNo);
213 
222 Vector InitAverage(real x, real y=0., real z=0.);
223 
232 real InitResistivity(real x, real y=0., real z=0.);
233 
240 void InitParameters();
241 
247 void InitTimeStep();
248 
256 void InitTree(Node* Root);
257 
265 Vector Limiter(const Vector u, const Vector v);
266 
273 real Limiter(const real x);
274 
282 real MinAbs(const real a, const real b);
283 
291 real NormMaxQuantities(const Vector& V);
292 
300 void Performance(const char* FileName);
301 
308 void PrintIntegral(const char* FileName);
309 
318 void RefreshTree(Node* Root);
319 
327 void Remesh(Node* Root);
328 
335 Vector FluxX(const Vector& Avg);
336 
343 Vector FluxY(const Vector& Avg);
344 
351 Vector FluxZ(const Vector& Avg);
352 
366 Vector SchemeHLL(const Cell& Cell1, const Cell& Cell2, const Cell& Cell3, const Cell& Cell4, const int AxisNo);
367 
381 Vector SchemeHLLD(const Cell& Cell1, const Cell& Cell2, const Cell& Cell3, const Cell& Cell4, const int AxisNo);
382 
398 Matrix stateUstar(const Vector& AvgL, const Vector& AvgR, const real prel, const real prer, real& slopeLeft, real& slopeRight,real& slopeM, real& slopeLeftStar, real& slopeRightStar, int AxisNo);
399 
409 void fluxCorrection(Vector& Flux,const Vector& AvgL, const Vector& AvgR, int AxisNo);
410 
418 void ShowTime(Timer arg);
419 
426 int Sign(const real a);
427 
436 Vector ArtificialViscosity(const Vector& Cell1, const Vector& Cell2, real dx, int AxisNo);
437 
448 Vector ResistiveTerms(Cell& Cell1, Cell& Cell2, Cell& Cell3, Cell& Cell4, int AxisNo);
449 
456 Vector Source(Cell& UserCell);
457 
464 real Step(real x);
465 
473 void TimeEvolution(FineMesh* Root);
474 
482 void TimeEvolution(Node* Root);
483 
493 void View(FineMesh* Root, const char* AverageFileName);
494 
506 void View(Node* Root, const char* TreeFileName, const char* MeshFileName, const char* AverageFileName);
507 
516 void ViewEvery(FineMesh* Root, int arg);
517 
528 void ViewEvery(Node* Root, int arg);
529 
537 void ViewIteration(FineMesh* Root);
538 
548 void ViewIteration(Node* Root);
549 
550 
551 //Parallel
552 void CreateMPIType(FineMesh *Root);
553 
561 void CPUExchange(FineMesh *Root,int);
562 
569 void FreeMPIType();
570 
577 void CreateMPITopology();
578 
585 void CreateMPILinks();
586 
593 void ReduceIntegralValues();
594 
595 
596 
597 
void Backup(Node *Root)
Stores the tree structure and data in order to restart a multiresolution computation.
Definition: Backup.cpp:30
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
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
real InitResistivity(real x, real y=0., real z=0.)
Returns the initial resistivity condition in (x, y, z) form the one defined in carmen.eta.
Definition: InitAverage.cpp:50
An object Node is an element of a graded tree structure, used for multiresolution computations...
Definition: Node.h:38
void FreeMPIType()
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:247
int DigitNumber(int arg)
Returns the number of digits of the integer arg.
Definition: DigitNumber.cpp:22
Vector FluxY(const Vector &Avg)
Returns the physical flux of MHD equations in Y direction.
Definition: PhysicalFluxMHD.cpp:58
An object FineMesh is a regular fine mesh, used for finite volume computations. It is not used for mu...
Definition: FineMesh.h:40
Vector ArtificialViscosity(const Vector &Cell1, const Vector &Cell2, real dx, int AxisNo)
Returns the artificial diffusion source terms in the cell UserCell.
Definition: ArtificialViscosity.cpp:11
Standard class for every vector in Carmen.
Definition: Vector.h:29
void Remesh(Node *Root)
Remesh the tree structure after a time evolution. The root node is Root. Only for multiresolution com...
Definition: Remesh.cpp:22
void InitTree(Node *Root)
Inits tree structure from initial condition, starting form the node Root. Only for multiresolution co...
Definition: InitTree.cpp:22
void PrintIntegral(const char *FileName)
Writes the integral values, like e.g flame velocity, global error, into file FileName.
Definition: PrintIntegral.cpp:30
int BoundaryRegion(const Vector &X)
Returns the boundary region type at the position X=(x,y,z). The returned value correspond to: 0 = Fl...
Definition: BoundaryRegion.cpp:63
int BC(int i, int AxisNo, int N=(1<< ScaleNb))
Returns the position of i taking into account the boundary conditions in the direction AxisNo...
Definition: BC.cpp:37
int FileWrite(FILE *f, const char *format, real arg)
Writes in binary or ASCII mode the real number arg into the file f with the format format...
Definition: FileWrite.cpp:22
void CreateMPIType(FineMesh *Root)
Definition: Parallel.cpp:121
void GetBoundaryCells(Cell &Cell1, Cell &Cell2, Cell &Cell3, Cell &Cell4, Cell &C1, Cell &C2, Cell &C3, Cell &C4, int AxisNo)
Transform the 4 cells of the flux Cell1, Cell2, Cell3, Cell4 into C1, C2, C3, C4, to take into accoun...
Definition: GetBoundaryCells.cpp:24
Vector Source(Cell &UserCell)
Returns the source term in the cell UserCell.
Definition: Source.cpp:23
real Step(real x)
Returns a step (1 if x <0, 0 if x >0, 0.5 if x=0)
Definition: Step.cpp:24
void TimeEvolution(FineMesh *Root)
Computes a time evolution on the regular fine mesh Root. Only for finite volume computations.
Definition: TimeEvolution.cpp:77
Vector Limiter(const Vector u, const Vector v)
Returns the value of the slope limiter between the slopes u and v.
Definition: Limiter.cpp:56
Vector FluxX(const Vector &Avg)
Returns the physical flux of MHD equations in X direction.
Definition: PhysicalFluxMHD.cpp:8
void CreateMPILinks()
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:271
void ViewIteration(FineMesh *Root)
Same as previous for a fine mesh Root. Only for finite volume computations.
Definition: ViewIteration.cpp:60
real MinAbs(const real a, const real b)
Returns the minimum in module of a and b.
Definition: MinAbs.cpp:22
void ViewEvery(FineMesh *Root, int arg)
Same as previous for a fine mesh Root. Only for finite volume computations.
Definition: ViewEvery.cpp:53
void AdaptTimeStep()
Adapts time step when required.
Definition: AdaptTimeStep.cpp:24
void fluxCorrection(Vector &Flux, const Vector &AvgL, const Vector &AvgR, int AxisNo)
This function apply the divergence-free correction to the numerical flux.
Definition: FluxCorrection.cpp:9
int Sign(const real a)
Returns 1 if a is non-negative, -1 elsewhere.
Definition: Sign.cpp:22
void CreateMPITopology()
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:22
Standard class for every matrix in Carmen.
Definition: Matrix.h:28
Vector FluxZ(const Vector &Avg)
Returns the physical flux of MHD equations in Z direction.
Definition: PhysicalFluxMHD.cpp:106
Vector ResistiveTerms(Cell &Cell1, Cell &Cell2, Cell &Cell3, Cell &Cell4, int AxisNo)
Returns the resistive source terms in the cell UserCell.
Definition: ResistiveTerms.cpp:11
Vector Flux(Cell &Cell1, Cell &Cell2, Cell &Cell3, Cell &Cell4, int AxisNo)
Returns the flux at the interface between Cell2 and Cell3. Here a 4-point space scheme is used...
Definition: Flux.cpp:22
void ReduceIntegralValues()
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:469
void ShowTime(Timer arg)
Writes on screen the estimation of total and remaining CPU times. These informations are stored in th...
Definition: ShowTime.cpp:23
void Performance(const char *FileName)
Computes the performance of the computation and, for cluster computations, write it into file FileNam...
Definition: Performance.cpp:22
Vector SchemeHLL(const Cell &Cell1, const Cell &Cell2, const Cell &Cell3, const Cell &Cell4, const int AxisNo)
Returns the HLL numerical flux for MHD equations. The scheme uses four cells to estimate the flux at ...
Definition: SchemeHLL.cpp:11
void InitTimeStep()
Inits time step and all the parameters which depend on it.
Definition: InitTimeStep.cpp:22
double CPUTimeRef(int iterations, int scales)
Returns the time required by a finite volume computation using iterations iterations and scales scale...
Definition: CPUTimeRef.cpp:23
real ComputedTolerance(const int ScaleNo)
Returns the computed tolerance at the scale ScaleNo, either using Harten or Donoho thresholding (if C...
Definition: ComputedTolerance.cpp:22
void View(FineMesh *Root, const char *AverageFileName)
Writes the current cel–averages of the fine mesh Root into file AverageFileName. Only for finite volu...
Definition: View.cpp:87
real NormMaxQuantities(const Vector &V)
Returns the Max-norm of the vector where every quantity is divided by its characteristic value...
Definition: NormMaxQuantities.cpp:25
void InitParameters()
Inits parameters from file carmen.par. If a parameter is not mentioned in this file, the default value is used.
Definition: Parameters.cpp:283
An object Timer gives information on the CPU time of long-time computations.
Definition: Timer.h:31
Vector SchemeHLLD(const Cell &Cell1, const Cell &Cell2, const Cell &Cell3, const Cell &Cell4, const int AxisNo)
Returns the HLLD numerical flux for MHD equations. The scheme uses four cells to estimate the flux at...
Definition: SchemeHLLD.cpp:12
Vector InitAverage(real x, real y=0., real z=0.)
Returns the initial condition in (x, y, z) form the one defined in carmen.ini.
Definition: InitAverage.cpp:22
void CPUExchange(FineMesh *Root, int)
Parallel function DOES NOT WORK!
Definition: Parallel.cpp:350
Matrix stateUstar(const Vector &AvgL, const Vector &AvgR, const real prel, const real prer, real &slopeLeft, real &slopeRight, real &slopeM, real &slopeLeftStar, real &slopeRightStar, int AxisNo)
Returns the intermediary states of HLLD numerical flux for MHD equations.
Definition: IntermediaryStates.cpp:12
#define real
Definition: PreProcessor.h:31
This header contains all parameters as global variables.