Peano
Loading...
Searching...
No Matches
CCZ4KernelTest.cpp
Go to the documentation of this file.
1#include "CCZ4KernelTest.h"
2
4#include "peano4/utils/Loop.h"
5#include "tarch/la/Vector.h"
6
7#include "tarch/logging/Log.h"
8
10//#include "exahype2/CellData.h"
11//#include "exahype2/fd/Functors.h"
12#include "exahype2/fd/fd4/FD4.h"
13#include "../applications/exahype2/ccz4/CCZ4Kernels.h"
14#include "../applications/exahype2/ccz4/CCZ4Kernels.cpph"
15
19
20#include <algorithm>
21#include <iomanip>
22#include <math.h>
23
24
26 TestCase ("exahype2::fd::tests::CCZ4KernelTest") {
27}
28
29
31 testMethod (AppleWithAppleTest);
32}
33
34
36 //in this test, the out layers inside the patch is not valid as the volumes
37 // in the halos can not update their first-derivatives. So only volumes with index between [3,26]^3
38 // holds valid test data.
39
40 // make no sense to test in 2d
41 #if Dimensions==3
42 constexpr int unknowns = 59;
43 constexpr int auxiliaryVariables = 0;
44 constexpr int numberOfDoFsPerAxisInCell = 30;
45 constexpr int haloSize = 3;
46 constexpr int numberOfDoFsPerAxisInCellWithHalos= 30+2*haloSize;
47
48 const tarch::la::Vector<3, double> x({0.5, 0.5, 0.5});
49 const tarch::la::Vector<3, double> cellSize({30.0/29.0, 30.0/29.0, 30.0/29.0});
50 tarch::la::Vector<3, int> currentPosition({0,0,0});
51 tarch::la::Vector<3, double> volumeCoordinates({0,0,0});
52
53 double oldQWithHalos[numberOfDoFsPerAxisInCellWithHalos*numberOfDoFsPerAxisInCellWithHalos*numberOfDoFsPerAxisInCellWithHalos*(unknowns+auxiliaryVariables)]={0};
54 double QOut[numberOfDoFsPerAxisInCell*numberOfDoFsPerAxisInCell*numberOfDoFsPerAxisInCell*(unknowns)]={0};
55
56 //random assignment here, keep the same with GRChombo.
57 const int scenario=1; // 0-flat spacetime test case, 1-apple with apple test
58
59 if (scenario==0){ //flat space
60 for (int xx = -haloSize; xx < numberOfDoFsPerAxisInCell+haloSize; xx++)
61 for (int yy = -haloSize; yy < numberOfDoFsPerAxisInCell+haloSize; yy++)
62 for (int zz = -haloSize; zz < numberOfDoFsPerAxisInCell+haloSize; zz++){
63 currentPosition=gridCellIndex3d(xx,yy,zz); //current index, notice the halo is included
64 for (int i=0;i<3;i++){
65 volumeCoordinates(i) = x(i)-0.5*cellSize(i);
66 volumeCoordinates(i) += (double(currentPosition(i))+0.5)*cellSize(i)/numberOfDoFsPerAxisInCell; //current position
67 }
68 //notice if you index using this routie the index need to start from 0 not -halo!
69 int CellIndexSerialised = peano4::utils::dLinearised((currentPosition+tarch::la::Vector<Dimensions,int>(3)),
70 numberOfDoFsPerAxisInCellWithHalos);
71
72 //assign quantities
73 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 0] = 1.0; //\tilde{\gamma}_11
74 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 3] = 1.0; //\tilde{\gamma}_22
75 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 5] = 1.0; //\tilde{\gamma}_33
76 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+16] = 1.0; //\alpha
77 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+54] = 1.0; //\phi
78 }
79 } else if (scenario==1){ //apple with apple
80 for (int xx = -haloSize; xx < numberOfDoFsPerAxisInCell+haloSize; xx++)
81 for (int yy = -haloSize; yy < numberOfDoFsPerAxisInCell+haloSize; yy++)
82 for (int zz = -haloSize; zz < numberOfDoFsPerAxisInCell+haloSize; zz++){
83 currentPosition=gridCellIndex3d(xx,yy,zz); //current index, notice the halo is included
84 for (int i=0;i<3;i++){
85 volumeCoordinates(i) = x(i)-0.5*cellSize(i);
86 volumeCoordinates(i) += (double(currentPosition(i))+0.5)*cellSize(i)/numberOfDoFsPerAxisInCell; //current position
87 }
88
89 int CellIndexSerialised = peano4::utils::dLinearised((currentPosition+tarch::la::Vector<Dimensions,int>(3)),
90 numberOfDoFsPerAxisInCellWithHalos);
91 //assign quantities
92 double **g, **K;
93 double phi, trK;
94 double Theta, Gamma[3];
95 double alpha, beta[3], b[3];
96 g = new double* [3]; K = new double* [3];
97 for (int l=0;l<3;l++){
98 g[l]=new double[3]; K[l]=new double[3];
99 }
100
101
102 prepareFieldData(g, phi, K, trK, Theta, Gamma, alpha, beta, b,
103 volumeCoordinates(0), volumeCoordinates(1), volumeCoordinates(2));
104
105 //for (int i=0;i<3;i++)
106 //for (int j=0;j<3;j++)
107 //double a=1.0/NAN;
108 //std::cout<<1.0/NAN<<std::endl;
109
110 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 0] = g[0][0]; // tilde{gamma}_ij
111 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 1] = g[0][1];
112 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 2] = g[0][2];
113 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 3] = g[1][1];
114 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 4] = g[1][2];
115 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 5] = g[2][2];
116
117 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 6] = K[0][0]; // tilde{A}_ij
118 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 7] = K[0][1];
119 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 8] = K[0][2];
120 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+ 9] = K[1][1];
121 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+10] = K[1][2];
122 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+11] = K[2][2];
123
124 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+12] = Theta; // Theta
125 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+13] = Gamma[0]; // hat{Gamma}^i
126 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+14] = Gamma[1];
127 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+15] = Gamma[2];
128
129 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+16] = alpha; // alpha
130 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+17] = beta[0]; // beta^i
131 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+18] = beta[1];
132 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+19] = beta[2];
133 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+20] = b[0]; // b^i
134 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+21] = b[1];
135 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+22] = b[2];
136
137 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+53] = trK; // K
138 oldQWithHalos[CellIndexSerialised*(unknowns+auxiliaryVariables)+54] = phi; // phi
139
140 delete g;
141 delete K;
142 }
143 }
144
145 ::exahype2::CellData patchData(oldQWithHalos, x, cellSize, 0, 0.01, QOut);
146
147 //set parameters
148 static constexpr double CCZ4KOSigma=8.0;
149 static constexpr double CCZ4itau = 1.0;
150 static constexpr double CCZ4k1 = 0.1;
151 static constexpr double CCZ4k2 = 0.0;
152 static constexpr double CCZ4k3 = 0.5;
153 static constexpr double CCZ4eta = 1.0;
154 static constexpr double CCZ4f = 0.75;
155 static constexpr double CCZ4g = 0.0;
156 static constexpr double CCZ4xi = 1.0;
157 static constexpr double CCZ4e = 1.0;
158 static constexpr double CCZ4c = 1.0;
159 static constexpr double CCZ4mu = 0.2;
160 static constexpr double CCZ4ds = 1.0;
161 static constexpr double CCZ4sk = 1.0;
162 static constexpr double CCZ4bs = 1.0;
163 static constexpr int CCZ4LapseType = 1;
164 static constexpr int CCZ4SO = 1;
165
166 //then calculate the derivatives first
167 //there is a potential seg error if number of auxiliary variables is not zero.
168 //but without auxiliary variables everything should be fine.
169 //notice only the volumes inside the patch get updated
171 patchData,
172 numberOfDoFsPerAxisInCell,
173 haloSize,
174 unknowns,
175 auxiliaryVariables
176 );
177
178 //call the kernel
180 patchData,
181 numberOfDoFsPerAxisInCell,
182 haloSize,
183 unknowns,
184 auxiliaryVariables,
185 CCZ4KOSigma,
186 false ,
187 true ,
188 true ,
189 false, // do not copy the old solution. Abuse this to get the RHS directly
191 [&](
192 const double * __restrict__ Q,
193 const tarch::la::Vector<Dimensions,double>& faceCentre,
195 double t,
196 double dt,
197 int normal,
198 double * __restrict__ F
199 )->void {
200
201 },
202 [&](
203 const double * __restrict__ Q,
204 const double * __restrict__ deltaQ,
205 const tarch::la::Vector<Dimensions,double>& faceCentre,
207 double t,
208 double dt,
209 int normal,
210 double * __restrict__ BTimesDeltaQ
211 )->void {
212 double gradQSerialised[unknowns*3];
213 for (int i=0; i<unknowns; i++) {
214 gradQSerialised[i+0*unknowns] = 0.0;
215 gradQSerialised[i+1*unknowns] = 0.0;
216 gradQSerialised[i+2*unknowns] = 0.0;
217
218 gradQSerialised[i+normal*unknowns] = deltaQ[i];
219 }
220 //for (int i=NumberOfUnknowns; i<NumberOfUnknowns+NumberOfAuxiliaryVariables; i++) {
221 // BgradQ[i] = 0.0;
222 //}
223 applications::exahype2::ccz4::ncp(BTimesDeltaQ, Q, gradQSerialised, normal, CCZ4LapseType, CCZ4ds, CCZ4c, CCZ4e, CCZ4f, CCZ4bs, CCZ4sk, CCZ4xi, CCZ4mu, CCZ4SO);
224 },
225 [&](
226 const double * __restrict__ Q,
229 double t,
230 double dt,
231 double * __restrict__ S
232 )->void {
233 tarch::memset(S, 0, unknowns*sizeof(double));
234 applications::exahype2::ccz4::source(S,Q, CCZ4LapseType, CCZ4ds, CCZ4c, CCZ4e, CCZ4f, CCZ4bs, CCZ4sk, CCZ4xi, CCZ4itau, CCZ4eta, CCZ4k1, CCZ4k2, CCZ4k3, CCZ4SO);
235 }
236 );
237
238 //examine the output, notice voulme index should between [3,26]^3
239 //tarch::la::Vector<3, int> volumeIndex({15,15,15});
240 //tarch::la::Vector<3, int> volumeIndex({8,20,17});
242
243 //15,15,15
244 //double dataArrayOfGRChombo[25]={-1.01737, -0.952962, 1.5358, 2.07569, -0.274332, -1.34847, 4.52609, 1.69924, 0.581672, -1.01265, 1.01734, -3.28122, -2.00355, 3.42196, 23.2832, 0.698199, 26.806, -2.87506, -1.21849, -4.32589, 4.44389, 18.9259, -2.17062, -11.1472, -0.0253941};
245 //8,20,17
246 //double dataArrayOfGRChombo[25]={-1.39787, -0.754864, 1.71128, 1.41365, 1.63204, -0.675397, 6.08035, 2.31647, 1.3964, -1.59352, 1.60353, -4.44641, -0.945263, 2.65077, 36.1263, 2.48858, 36.1493, -5.74057, -4.27386, -8.90803, 3.83137, 29.8071, -3.80586, -10.7919, -0.0296434};
247 //23,9,18
248 double dataArrayOfGRChombo[25]={-3.1866, 0.365012, 3.35262, 0.0372327, 3.36914, 1.34624, 4.88302, 3.09191, 0.296349, -3.06936, 0.936401, -2.12479, 3.00939, 9.47752, 35.5824, -9.94656, 42.4051, -0.753163, 3.48535, -2.7644, 10.7661, 34.7114, -11.501, -9.21531, -0.359624};
249
250 int checkingQ=0;
251 exahype2::enumerator::AoSLexicographicEnumerator QInEnumerator (1,numberOfDoFsPerAxisInCell,haloSize,unknowns,auxiliaryVariables);
252 exahype2::enumerator::AoSLexicographicEnumerator QOutEnumerator(1,numberOfDoFsPerAxisInCell,0,unknowns,0); //halosize and auxiliary variables are both 0,
253 for (int i=0;i<3;i++){
254 volumeCoordinates(i) = patchData.cellCentre[0](i)-0.5*patchData.cellSize[0](i);
255 volumeCoordinates(i) += (double(volumeIndex(i))+0.5)*patchData.cellSize[0](i)/QOutEnumerator._numberOfDoFsPerAxisInCell;
256 }
257
258 /*validateNumericalEqualsWithParams3(0, 1,
259 patchData.toString(),
260 patchData.QOut[0][QOutEnumerator(0, volumeIndex, checkingQ)],
261 volumeCoordinates
262 );*/
263 //for (checkingQ=0;checkingQ<unknowns;checkingQ++){
264 // std::cout<<patchData.QOut[0][QOutEnumerator(0, volumeIndex, checkingQ)]<<"\t";
265 //}
266 //std::cout<<std::endl;
267 //std::cout << "Index: " << volumeIndex << " Coordinates: " << volumeCoordinates << std::endl;
268 //std::cout << std::fixed <<std::setprecision(6);
269 /*for (checkingQ=0;checkingQ<23;checkingQ++){
270 std::cout<<"dtQ["<<checkingQ<<"]\t= "<<patchData.QOut[0][QOutEnumerator(0, volumeIndex, checkingQ)]<<std::endl;
271 }
272 std::cout<<"dtQ[53]\t= "<<patchData.QOut[0][QOutEnumerator(0, volumeIndex, 53)]<<std::endl;
273 std::cout<<"dtQ[54]\t= "<<patchData.QOut[0][QOutEnumerator(0, volumeIndex, 54)]<<std::endl;
274 std::cout<<"dtchi=2*Q[54]*dtQ[54]\t= "
275 <<2*patchData.QIn[0][QInEnumerator(0, volumeIndex, 54)]*patchData.QOut[0][QOutEnumerator(0, volumeIndex, 54)]
276 <<std::endl;*/
277
278 /*for (checkingQ=0;checkingQ<23;checkingQ++){
279 std::cout<<"dtQ["<<checkingQ<<"]\t= "<<patchData.QOut[0][QOutEnumerator(0, volumeIndex, checkingQ)]
280 <<"\tGRChombo result = "<<dataArrayOfGRChombo[checkingQ]<<"\t Ratio = "
281 <<100*patchData.QOut[0][QOutEnumerator(0, volumeIndex, checkingQ)]/dataArrayOfGRChombo[checkingQ]
282 <<"%"<<std::endl;
283 }
284
285 std::cout<<"dtQ[53]\t= "<<patchData.QOut[0][QOutEnumerator(0, volumeIndex, 53)]
286 <<"\tGRChombo result = "<<dataArrayOfGRChombo[23]<<"\t Ratio = "
287 <<100*patchData.QOut[0][QOutEnumerator(0, volumeIndex, 53)]/dataArrayOfGRChombo[23]
288 <<"%"<<std::endl;
289
290 double dtphi=patchData.QOut[0][QOutEnumerator(0, volumeIndex, 54)];
291 std::cout<<"dtchi=2*Q[54]*dtQ[54]\t= "
292 <<2*patchData.QIn[0][QInEnumerator(0, volumeIndex, 54)]*dtphi
293 <<"\tGRChombo result = "<<dataArrayOfGRChombo[24]<<"\t Ratio = "
294 <<100*2*patchData.QIn[0][QInEnumerator(0, volumeIndex, 54)]*dtphi/dataArrayOfGRChombo[24]
295 <<"%"<<std::endl;*/
296
297 /*for (checkingQ=0;checkingQ<23;checkingQ++){
298 validateNumericalEqualsWithParams3(0, 1,
299 patchData.toString(),
300 patchData.QOut[0][QOutEnumerator(0, volumeIndex, checkingQ)],
301 volumeCoordinates
302 );
303 }*/
304
305
306 #endif
307}
308
309
311 double** g,
312 double& phi,
313 double** K,
314 double& trK,
315 double& Theta,
316 double* Gamma,
317 double& alpha,
318 double* beta,
319 double* b,
320 double x, double y, double z
321){
322 double g_UU[3][3];
323
324 g[0][0] = 1.36778 + 2.39731 * x + 4.53541 * x * x +
325 19.9771 * x * y * y * y + 6.13801 * y * z +
326 5.65185 * z * z + 9.35842 * z * z * z * z;
327 g[0][1] = -0.07646 - 0.48786 * x - 0.75098 * x * x -
328 1.73683 * x * y * y * y + 1.71676 * y * z +
329 1.03662 * z * z + 0.35630 * z * z * z * z;
330 g[0][2] = -0.10083 + 0.12445 * x - 1.26649 * x * x -
331 1.95052 * x * y * y * y + 0.73091 * y * z -
332 1.49835 * z * z - 2.39024 * z * z * z * z;
333 g[1][1] = 0.84072 + 2.31163 * x + 3.32275 * x * x +
334 15.1662 * x * y * y * y + 8.48730 * y * z +
335 3.05098 * z * z + 17.8448 * z * z * z * z;
336 g[1][2] = -0.42495 - 0.33464 * x - 0.47012 * x * x -
337 7.38477 * x * y * y * y + 0.41896 * y * z -
338 1.36394 * z * z + 5.25894 * z * z * z * z;
339 g[2][2] = 0.60995 + 1.30428 * x + 3.86237 * x * x +
340 22.7614 * x * y * y * y + 6.93818 * y * z +
341 4.39250 * z * z + 19.0244 * z * z * z * z;
342 g[1][0] = g[0][1];
343 g[2][0] = g[0][2];
344 g[2][1] = g[1][2];
345
346 double detg = g[0][0] * g[1][1] * g[2][2] +
347 2 * g[0][1] * g[0][2] * g[1][2] -
348 g[0][0] * g[1][2] * g[1][2] -
349 g[1][1] * g[0][2] * g[0][2] -
350 g[2][2] * g[0][1] * g[0][1];
351 g_UU[0][0] = (g[1][1] * g[2][2] - g[1][2] * g[1][2]) / detg;
352 g_UU[0][1] = (g[0][2] * g[1][2] - g[0][1] * g[2][2]) / detg;
353 g_UU[0][2] = (g[0][1] * g[1][2] - g[0][2] * g[1][1]) / detg;
354 g_UU[1][1] = (g[0][0] * g[2][2] - g[0][2] * g[0][2]) / detg;
355 g_UU[1][2] = (g[0][1] * g[0][2] - g[0][0] * g[1][2]) / detg;
356 g_UU[2][2] = (g[0][0] * g[1][1] - g[0][1] * g[0][1]) / detg;
357 g_UU[1][0] = g_UU[0][1];
358 g_UU[2][0] = g_UU[0][2];
359 g_UU[2][1] = g_UU[1][2];
360
361 phi=std::pow(detg,-1./6.);
362 double phisq=phi*phi;
363 trK=0;
364 K[0][0] = -0.16238 - 0.74295 * x + 0.51595 * x * x -
365 6.60239 * x * y * y * y - 0.76401 * y * z -
366 1.81131 * z * z - 3.88228 * z * z * z * z;
367 K[0][1] = 0.15054 - 0.60088 * x - 0.15428 * x * x +
368 3.16779 * x * y * y * y - 2.00687 * y * z -
369 1.35442 * z * z - 0.67601 * z * z * z * z;
370 K[0][2] = -0.02174 - 0.36243 * x + 0.81531 * x * x +
371 4.34918 * x * y * y * y + 0.90419 * y * z -
372 0.85088 * z * z - 6.45097 * z * z * z * z;
373 K[1][1] = -0.47653 - 0.43889 * x + 0.87342 * x * x +
374 4.24684 * x * y * y * y + 0.26290 * y * z +
375 1.90095 * z * z + 3.69515 * z * z * z * z;
376 K[1][2] = 0.37472 + 0.03657 * x - 0.10327 * x * x -
377 0.95744 * x * y * y * y - 1.20800 * y * z -
378 0.43064 * z * z - 0.25419 * z * z * z * z;
379 K[2][2] = 0.34184 + 0.21495 * x - 0.73195 * x * x +
380 7.81626 * x * y * y * y + 2.48359 * y * z +
381 1.89657 * z * z - 4.10980 * z * z * z * z;
382 K[1][0] = K[0][1];
383 K[2][0] = K[0][2];
384 K[2][1] = K[1][2];
385 for (int i=0;i<3;i++)
386 for (int j=0;j<3;j++) trK += g_UU[i][j]*K[i][j];
387
388 for (int i=0;i<3;i++)
389 for (int j=0;j<3;j++) {
390 g[i][j]=phisq*g[i][j]; //\tilde{\gamma}_ij
391 K[i][j]=phisq*(K[i][j]- 1./3. * trK * g[i][j]/phisq); //\tilde{A}_ij
392 }
393
394 Theta = 0.27579 + 0.25791 * x + 1.40488 * x * x +
395 5.68276 * x * y * y * y +
396 3.04325 * y * z + 1.81250 * z * z +
397 1.01832 * z * z * z * z;
398 Gamma[0] = -0.49482 + 0.89227 * x + 0.05571 * x * x -
399 5.38570 * x * y * y * y + 0.13979 * y * z -
400 0.68588 * z * z - 4.39964 * z * z * z * z;
401 Gamma[1] = -0.09082 - 0.31017 * x + 1.06980 * x * x +
402 7.81524 * x * y * y * y - 1.65016 * y * z -
403 0.53352 * z * z - 3.20997 * z * z * z * z;
404 Gamma[2] = -0.42367 + 0.03891 * x - 0.87898 * x * x +
405 6.67657 * x * y * y * y - 3.44662 * y * z -
406 0.19655 * z * z + 2.97524 * z * z * z * z;
407
408 alpha = 0.73578 + 0.36898 * x + 0.64348 * x * x +
409 9.33487 * x * y * y * y +
410 0.99469 * y * z + 0.20515 * z * z +
411 8.88385 * z * z * z * z;
412
413
414 beta[0] = 0.00000 + 0.18795 * x - 0.52389 * x * x -
415 4.14079 * x * y * y * y +
416 0.73135 * y * z - 0.27057 * z * z +
417 3.24187 * z * z * z * z;
418 beta[1] = 0.00000 - 0.30316 * x - 0.15184 * x * x -
419 0.48815 * x * y * y * y +
420 2.45991 * y * z - 0.79248 * z * z +
421 7.14007 * z * z * z * z;
422 beta[2] = 0.00000 + 0.68835 * x - 0.52219 * x * x -
423 7.50449 * x * y * y * y -
424 2.35372 * y * z - 0.21476 * z * z +
425 4.36363 * z * z * z * z;
426
427 b[0] = -0.26928 + 0.35045 * x - 0.48884 * x * x +
428 2.72465 * x * y * y * y - 2.59022 * y * z -
429 0.27384 * z * z + 0.38748 * z * z * z * z;
430 b[1] = 0.40234 + 0.26741 * x + 1.94822 * x * x -
431 0.78276 * x * y * y * y + 2.12346 * y * z +
432 0.69086 * z * z - 4.47639 * z * z * z * z;
433 b[2] = 0.40313 + 0.00569 * x - 1.12452 * x * x -
434 5.49255 * x * y * y * y - 2.21932 * y * z +
435 0.49523 * z * z + 1.29460 * z * z * z * z;
436}
#define testMethod(name)
Run a test method and check for errors.
Definition TestMacros.h:24
virtual void run()
This routine is triggered by the TestCaseCollection.
void prepareFieldData(double **g, double &phi, double **K, double &trK, double &Theta, double *Gamma, double &alpha, double *beta, double *b, double x, double y, double z)
KeywordToAvoidDuplicateSymbolsForInlinedFunctions void source(double *S, const double *const Q, const int CCZ4LapseType, const double CCZ4ds, const double CCZ4c, const double CCZ4e, const double CCZ4f, const double CCZ4bs, const double CCZ4sk, const double CCZ4xi, const double CCZ4itau, const double CCZ4eta, const double CCZ4k1, const double CCZ4k2, const double CCZ4k3, const double CCZ4SO) InlineMethod
The source term is one out of two terms that we use in our CCZ4 formulation.
KeywordToAvoidDuplicateSymbolsForInlinedFunctions void ncp(double *BgradQ, const double *const Q, const double *const gradQSerialised, const int normal, const int CCZ4LapseType, const double CCZ4ds, const double CCZ4c, const double CCZ4e, const double CCZ4f, const double CCZ4bs, const double CCZ4sk, const double CCZ4xi, const double CCZ4mu, const double CCZ4SO) InlineMethod
void reconstruct_first_derivatives(::exahype2::CellData< double, double > &patchData, int numberOfGridCellsPerPatchPerAxis, int haloSize, int unknowns, int auxiliaryVariables)
Helper routine to reconstruct the first derivatives.
void timeStep_patchwise_heap_functors(::exahype2::CellData< double, double > &patchData, int numberOfGridCellsPerPatchPerAxis, int haloSize, int unknowns, int auxiliaryVariables, double KOSigma, bool evaluateFlux, bool evaluateDifferentialSource, bool evaluateAlgebraicSource, bool copyOldTimeStepAndScaleWithTimeStepSize, DifferentialSourceTermVariant variant, Flux flux, NonconservativeProduct DifferentialSource, Source AlgebraicSource)
Fourth-order Finite Differences.
auto volumeIndex(Args... args)
Definition VolumeIndex.h:54
CPUGPUMethod int dLinearised(const tarch::la::Vector< Dimensions, int > &counter, int max)
Map d-dimensional vector onto integer index.
Definition Loop.cpp:106
double * memset(double *dest, double ch, size_t byteCount)
Alternative GPU-ready version of memset.
const struct part *restrict const struct part *restrict const float const float beta
Representation of a number of cells which contains all information that's required to process the sto...
Definition CellData.h:77
tarch::la::Vector< Dimensions, double > * cellCentre
Definition CellData.h:83
tarch::la::Vector< Dimensions, double > * cellSize
Definition CellData.h:84
Simple vector class.
Definition Vector.h:150