-
Notifications
You must be signed in to change notification settings - Fork 9
/
timeintegration.h
134 lines (105 loc) · 3.36 KB
/
timeintegration.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
/**
* @author Christoph Schaefer [email protected] and Thomas I. Maindl
*
* @section LICENSE
* Copyright (c) 2019 Christoph Schaefer
*
* This file is part of miluphcuda.
*
* miluphcuda is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* miluphcuda is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with miluphcuda. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef _TIMEINTEGRATION_H
#define _TIMEINTEGRATION_H
#include <stdio.h>
#include <libconfig.h>
#include <pthread.h>
#include <assert.h>
#include "parameter.h"
#include "miluph.h"
#include "io.h"
#include "cuda_utils.h"
#include "kernel.h"
#include "extrema.h"
#include "sinking.h"
// Courant (CFL) number (note that our sml is defined up to the zero of the kernel, not half of it)
#define COURANT_FACT 0.4
// factor for limiting timestep based on local forces/acceleration
#define FORCES_FACT 0.2
extern int startTimestep;
extern int numberOfTimesteps;
extern double startTime;
extern double timePerStep;
extern double dt_host;
extern double currentTime;
extern double h5time;
extern __device__ int numParticles;
extern __device__ int numPointmasses;
extern int *relaxedPerBlock;
extern void (*integrator)();
void timeIntegration(void);
// the available integrators
void euler(void);
void predictor_corrector(void);
void predictor_corrector_euler(void);
void rk2Adaptive(void);
void rk4_nbodies(void); // do we need it here? do we say.. bye,bye?
void heun_rk4(void);
double calculate_angular_momentum(void);
void copyToHostAndWriteToFile(int timestep, int lastTimestep);
__device__ int childListIndex(int nodeIndex, int childNumber);
__global__ void detectVelocityRelaxation(int *relaxedPerBlock);
__device__ int stressIndex(int particleIndex, int row, int col);
void afterIntegrationStep(void);
__global__ void symmetrizeStress(void);
#define NUM_THREADS_512 512
#define NUM_THREADS_256 256
#define NUM_THREADS_128 128
#define NUM_THREADS_64 64
#define NUM_THREADS_1 1
#define NUM_THREADS_COMPUTATIONAL_DOMAIN 128
#define NUM_THREADS_BUILD_TREE 32
#define NUM_THREADS_TREEDEPTH 128
#define NUM_THREADS_TREECHANGE 128
#define NUM_THREADS_CALC_CENTER_OF_MASS 256
#define NUM_THREADS_SELFGRAVITY 128
#define NUM_THREADS_BOUNDARY_CONDITIONS 128
#define NUM_THREADS_NEIGHBOURSEARCH 256
#define NUM_THREADS_SYMMETRIZE_INTERACTIONS 256
#define NUM_THREADS_DENSITY 256
#define NUM_THREADS_PRESSURE 256
#define NUM_THREADS_PALPHA_POROSITY 256
#define NUM_THREADS_REDUCTION 256
#define NUM_THREADS_DETECTRELAX 256
#define NUM_THREADS_LIMITTIMESTEP 256
// RK2
#define NUM_THREADS_RK2_INTEGRATE_STEP 256
#define NUM_THREADS_ERRORCHECK 256
// RK4
#define NUM_THREADS_RK4_INTEGRATE_STEP 256
// EULER
#define NUM_THREADS_EULER_INTEGRATOR 256
// PREDICTOR-CORRECTOR
#define NUM_THREADS_PC_INTEGRATOR 256
#define EMPTY -1
#define LOCKED -2
#define MAXDEPTH 128
// Runge-Kutta constants
#define B21 0.5
#define B31 -1.0
#define B32 2.0
#define C1 1.0
#define C2 4.0
#define C3 1.0
#endif