-
Notifications
You must be signed in to change notification settings - Fork 7
/
voigt.h
2342 lines (1911 loc) · 68.7 KB
/
voigt.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
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
//constants for sub-Lorentzian chi
__constant__ float sLchi_c[17];
__constant__ int useSubLorentzian_c[1];
// *********************************************
// This function calculates the B factors for the
// Perrin and Hartmann 1989 sub-Lorentzian chi factors
//
//Author Simon Grimm
//August 2020
// **********************************************
__host__ void subLorentzianB(double T){
float alpha1 = sLchi_h[3];
float beta1 = sLchi_h[4];
float epsilon1 = sLchi_h[5];
float alpha2 = sLchi_h[6];
float beta2 = sLchi_h[7];
float epsilon2 = sLchi_h[8];
float alpha3 = sLchi_h[9];
float beta3 = sLchi_h[10];
float epsilon3 = sLchi_h[11];
float B1 = alpha1 + beta1 * exp(-epsilon1 * T);
float B2 = alpha2 + beta2 * exp(-epsilon2 * T);
float B3 = alpha3 + beta3 * exp(-epsilon3 * T);
sLchi_h[14] = B1;
sLchi_h[15] = B2;
sLchi_h[16] = B3;
cudaMemcpyToSymbol(sLchi_c, sLchi_h, 17 * sizeof(float), 0, cudaMemcpyHostToDevice);
}
__host__ void subLorentzianConstantCopy(int useSubLorentzian){
cudaMemcpyToSymbol(useSubLorentzian_c, &useSubLorentzian, sizeof(int), 0, cudaMemcpyHostToDevice);
}
// *********************************************
// This function calculates the chi factors for the
// Perrin and Hartmann 1989 sub-Lorentzian correction
//
//Author Simon Grimm
//August 2020
// **********************************************
__device__ float sLChi(float Dnu){
float chi;
float sigma1 = sLchi_c[0];
float sigma2 = sLchi_c[1];
float sigma3 = sLchi_c[2];
float B1 = sLchi_c[14];
float B2 = sLchi_c[15];
float B3 = sLchi_c[16];
if(Dnu < sigma1){
chi = 1.0f;
}
else if(Dnu < sigma2){
chi = exp(-B1 * (Dnu - sigma1));
}
else if(Dnu < sigma3){
chi = exp(-B1 * (sigma2 - sigma1) - B2 * (Dnu - sigma2));
}
else{
chi = exp(-B1 * (sigma2 - sigma1) - B2 * (sigma3 - sigma2) - B3 * (Dnu - sigma3));
}
return chi;
}
// *********************************************
//This function calculates the Series Sigma1. Sigma2 and Sigma3 (Equations 27, 28, and 29) from Alg 916
//The parameter def_TOL sets a tolerance where to truncate the series
//It returns the values for sigma1 sigma2 and sigma3
//
// This implementation has still problems with loss of accuracy
//
//
//Author Simon Grimm, Adapted from Zaghloul & Ali, Algorithm 916
//November 2014
// **********************************************
__device__ void Sigma(const double x, const double y, double &s1, double &s2, double &s3, const double a, const double ex2, const int id){
s1 = 0.0;
double sold1 = s1;
s2 = 0.0;
double sold2 = s2;
s3 = 0.0;
double sold3 = s3;
double f, f3p, f3n;
double an, an3p, an3n;
double yy = y * y;
int n0 = (int)(ceil(x / a)); //starting point for sigma3 series
int n3p, n3n;
int stop1 = 0;
int stop2 = 0;
int stop3 = 0;
double ean2;
double e2axn = exp(-2.0 * a * x);
double e2axns = 1.0;
double an0 = a * n0;
double e3p = exp(-2.0 * a * (an0 - x - a));
double e3n = exp(2.0 * a * (an0 - x));
double e3ps = exp(-(an0 * an0 - 2.0 * a * an0 - 2.0 * an0 * x + x * x + a * a + 2.0 * a * x));
double e3ns = exp(-(an0 * an0 - 2.0 * an0 * x + x * x));
double st;
for(int n = 1; n < 100; ++n){
n3p = n0 + n - 1;
n3n = n0 - n;
an = a * n;
an3p = a * n3p;
an3n = a * n3n;
ean2 = exp(-an * an);
e2axns *= e2axn;
e3ps *= e3p;
e3ns *= e3n;
f = 1.0 / (an * an + yy);
f3p = 1.0 / (an3p * an3p + yy);
f3n = 1.0 / (an3n * an3n + yy);
st = f * ean2 * ex2;
s1 += st;
s2 += st * e2axns;
s3 += f3p * ean2 * e3ps;
if(n3n >= 1) s3 += f3n * ean2 * e3ns;
if(fabs(s1 - sold1) < def_TOL) stop1 = 1;
if(fabs(s2 - sold2) < def_TOL) stop2 = 1;
if(fabs(s3 - sold3) < def_TOL) stop3 = 1;
if(stop1 == 1 && stop2 ==1 && stop3 == 1) break;
sold1 = s1;
sold2 = s2;
sold3 = s3;
//if(n >= 100-1) printf("Sigma Series did not converge\n");
}
}
// *********************************************
//This function calculates the Series Sigma1. Sigma2 and Sigma3 (Equations 15, 16, and 17) from Alg 916
//The parameter def_TOL sets a tolerance where to truncate the series
//It returns the values for sigma1 sigma2 and sigma3
//Author Simon Grimm, Adapted from Zaghloul & Ali, Algorithm 916
//November 2014
// **********************************************
__device__ void Sigmab(double x, const double y, double &s1, double &s2, double &s3, const double a, const double ex2, const int id){
s1 = 0.0;
double sold1 = s1;
s2 = 0.0;
double sold2 = s2;
s3 = 0.0;
double sold3 = s3;
double f, f3p, f3n;
double an, an3p, an3n;
double yy = y * y;
if(x < 0.0) x = -x;
int n0 = (int)(ceil(x / a)); //starting point for sigma3 series
int n3p, n3n;
int stop1 = 0;
int stop2 = 0;
int stop3 = 0;
double e2axn = exp(-2.0 * a * x);
for(int n = 1; n < 100; ++n){
n3p = n0 + n - 1;
n3n = n0 - n;
an = a * n;
an3p = a * n3p;
an3n = a * n3n;
f = 1.0 / (an * an + yy);
f3p = 1.0 / (an3p * an3p + yy);
f3n = 1.0 / (an3n * an3n + yy);
s1 += f * exp(-(an * an + x * x));
s2 += f * exp(-(an + x) * (an + x));
s3 += f3p * exp(-(an3p - x) * (an3p - x));
if(n3n >= 1) s3 += f3n * exp(-(an3n - x) * (an3n - x));
if(fabs(s1 - sold1) < def_TOL) stop1 = 1;
if(fabs(s2 - sold2) < def_TOL) stop2 = 1;
if(fabs(s3 - sold3) < def_TOL) stop3 = 1;
if(stop1 == 1 && stop2 ==1 && stop3 == 1) break;
sold1 = s1;
sold2 = s2;
sold3 = s3;
// if(n >= 100-1) printf("Sigma Series did not converge %d\n", id);
}
}
__device__ void Sigmabf(float x, const float y, float &s1, float &s2, float &s3, const float a, const float ex2, const int id){
s1 = 0.0f;
float sold1 = s1;
s2 = 0.0f;
float sold2 = s2;
s3 = 0.0f;
float sold3 = s3;
float f, f3p, f3n;
float an, an3p, an3n;
float yy = y * y;
x = fabsf(x);
int n0 = (int)(ceil(x / a)); //starting point for sigma3 series
int n3p, n3n;
int stop1 = 0;
int stop2 = 0;
int stop3 = 0;
float e2axn = expf(-2.0f * a * x);
for(int n = 1; n < 100; ++n){
n3p = n0 + n - 1;
n3n = n0 - n;
an = a * n;
an3p = a * n3p;
an3n = a * n3n;
f = 1.0f / (an * an + yy);
f3p = 1.0f / (an3p * an3p + yy);
f3n = 1.0f / (an3n * an3n + yy);
s1 += f * expf(-(an * an + x * x));
s2 += f * expf(-(an + x) * (an + x));
s3 += f3p * expf(-(an3p - x) * (an3p - x));
if(n3n >= 1) s3 += f3n * expf(-(an3n - x) * (an3n - x));
if(fabs(s1 - sold1) < def_TOLF) stop1 = 1;
if(fabs(s2 - sold2) < def_TOLF) stop2 = 1;
if(fabs(s3 - sold3) < def_TOLF) stop3 = 1;
if(stop1 == 1 && stop2 ==1 && stop3 == 1) break;
sold1 = s1;
sold2 = s2;
sold3 = s3;
// if(n >= 100-1) printf("Sigma Series did not converge %d\n", id);
}
}
// *************************************************
//This function calculates the Voigt profile V(x,y) as equation 13 from Zaghloul & Ali, Algorithm 916
//it calls the Sigma function
//The parameter def_TOL sets a tolerance where to truncate the series
//Author Simon Grimm, Adapted from Zaghloul & Ali, Algorithm 916
//November 2014
// **************************************************
__device__ double voigt_916(const double x, const double y, const double a, const int id){
double s1, s2, s3;
double ex2 = exp(-x * x);
//Compute Sigma Series
if(x != 0.0 && y != 0.0) Sigmab(x, y, s1, s2, s3, a, ex2, id);
double xy = x * y;
double a2ipi = 2.0 * a / M_PI;
double cos2xy = cos(2.0 * xy);
double sinxy = sin(xy);
double t1 = ex2 * erfcx(y) * cos2xy;
t1 += a2ipi * x * sinxy * ex2 * sinxy / xy;
t1 += a2ipi * y * (-cos2xy * s1 + 0.5 * (s2 + s3));
if(x == 0) t1 = erfcx(y);
if(y == 0) t1 = ex2;
//if(x*x + y*y > 1.0e18) t1 = y / (sqrt(M_PI) * (x * x + y * y));
return t1;
}
__device__ float voigt_916f(const float x, const float y, const float a, const int id){
float s1, s2, s3;
float ex2 = expf(-x * x);
//Compute Sigma Series
if(x != 0.0f && y != 0.0f) Sigmabf(x, y, s1, s2, s3, a, ex2, id);
float xy = x * y;
float a2ipi = 2.0f * a / M_PIf;
float cos2xy = cosf(2.0f * xy);
float sinxy = sinf(xy);
float t1 = ex2 * erfcxf(y) * cos2xy;
t1 += a2ipi * x * sinxy * ex2 * sinxy / xy;
t1 += a2ipi * y * (-cos2xy * s1 + 0.5f * (s2 + s3));
if(x == 0) t1 = erfcxf(y);
if(y == 0) t1 = ex2;
//if(x*x + y*y > 1.0e18) t1 = y / (sqrt(M_PI) * (x * x + y * y));
return t1;
}
// *************************************************
//This kernel calculates the integrated line strength, the Lorentz and the Doppler halfwidths
//
//Author Simon Grimm
//October 2016
// *************************************************
__global__ void S2_kernel(double *nu_d, double *S_d, double *A_d, double *vy_d, double *ialphaD_d, double *n_d, double *delta_d, double *EL_d, int *ID_d, const int NL, const double T, const double P, const int kk){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx + kk;
if(id < NL){
double nu = nu_d[id] + delta_d[id] * P;
nu_d[id] = nu;
double alphaL = vy_d[id];
double ialphaD = ialphaD_d[id] / nu;
double EL = EL_d[id];
double c = def_h * def_c / def_kB;
S_d[id] *= exp(-EL * c / T + EL * c / def_T0) * (1.0 - exp(-c * nu / T)) / (1.0 - exp(-c * nu / def_T0)) * ialphaD;
ialphaD_d[id] = ialphaD; //inverse Doppler halfwidth
alphaL *= (P / def_POatm) * pow(def_T0 / T, n_d[id]);
alphaL += (A_d[id] / (4.0 * M_PI * def_c)); //1/cm
vy_d[id] = alphaL * ialphaD;
ID_d[id] = id;
if(nu == 0.0){
S_d[id] = 0.0;
ialphaD_d[id] = 0.0;
vy_d[id] = 0.0;
}
//if(id < 100) printf("%d %g %g %g %g %g\n", id, nu_d[id], S_d[id], ialphaD_d[id], EL, exp(-c * nu / T));
}
}
__global__ void Sf_kernel(double *nu_d, double *S_d, double *A_d, double *vy_d, double *ialphaD_d, double *n_d, double *EL_d, int *ID_d, const int NL, const double c, const double T1, const double P, const int kk){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx + kk;
if(id < NL){
double nu = nu_d[id];
double ialphaD = ialphaD_d[id] / nu;
double EL = EL_d[id];
S_d[id] *= exp(-c * EL) * (1.0 - exp(-c * nu)) * ialphaD;
ialphaD_d[id] = ialphaD; //inverse Doppler halfwidth
double alphaL = vy_d[id];
alphaL *= (P / def_PObar) * pow(T1, n_d[id]);
alphaL += A_d[id]; //1/cm
vy_d[id] = alphaL * ialphaD;
ID_d[id] = id;
if(nu == 0){
S_d[id] = 0.0;
ialphaD_d[id] = 0.0;
vy_d[id] = 0.0;
}
//if(id < 10) printf("ialphaD %d %g %g %g\n", id, nu_d[id], ialphaD_d[id], vy_d[id]);
//if( id < 100) printf("S %d %g\n", id, S_d[id]);
//if(id > 156000 && id < 158000) printf("S %d %g\n", id, S_d[id]);
}
}
__global__ void SfSuper_kernel(double *nu_d, double *S_d, double *vy_d, double *ialphaD_d, double *n_d, int *ID_d, const int NL, const double T1, const double P, const int kk){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx + kk;
if(id < NL){
double nu = nu_d[id];
double ialphaD = ialphaD_d[id] / nu;
S_d[id] *= ialphaD;
ialphaD_d[id] = ialphaD; //inverse Doppler halfwidth
double alphaL = vy_d[id];
alphaL *= (P / def_PObar) * pow(T1, n_d[id]);
vy_d[id] = alphaL * ialphaD;
ID_d[id] = id;
if(nu == 0){
S_d[id] = 0.0;
ialphaD_d[id] = 0.0;
vy_d[id] = 0.0;
}
//if(id < 10) printf("ialphaD %d %g %g %g\n", id, nu_d[id], ialphaD_d[id], vy_d[id]);
//if( id < 100) printf("S %d %g\n", id, S_d[id]);
//if(id > 156000 && id < 158000) printf("S %d %g\n", id, S_d[id]);
}
}
// *************************************************
//This kernel creates float arrays from double arrays and arrays which depends only on the pre-sorted other Line arrays
//It is called after the sorting routine to avoid sorting twice the same data
//Author: Simon Grimm
//Date: Oktober 2019
// *************************************************
__global__ void S3_kernel(double *nu_d, double *S_d, double *S1_d, double *vy_d, double *ialphaD_d, float* Sf_d, float *S1f_d, float *vyf_d, float *vcut2_d, float *va_d, float *vb_d, const double cut, const int cutMode, const int profile, const double numin, const double dnu, int useIndividualX, const int NL, const int kk){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx + kk;
if(id < NL){
double ialphaD = ialphaD_d[id];
vcut2_d[id] = (float)(cut * cut * ialphaD * ialphaD); //square of modified cut lenght
if(cutMode == 1){
vcut2_d[id] = (float)(cut * cut * vy_d[id] * vy_d[id]);
}
if(cutMode == 2){
vcut2_d[id] = (float)(cut * cut);
}
if(profile < 3){
S1_d[id] = S_d[id] * vy_d[id] / M_PI;
}
else if(profile == 3){
//Doppler
S1_d[id] = S_d[id] / sqrt(M_PI);
}
else{
//Binned Gaussian
S1_d[id] = S_d[id];
}
double nu = nu_d[id];
if(useIndividualX == 0){
va_d[id] = (float)((numin - nu) * ialphaD);
vb_d[id] = (float)(dnu * ialphaD);
//if(id == 304040) printf("vab %d %.20g %.20g %.20g %.20g\n", id, numin, nu, ialphaD, dnu);
//if(id == 304040) printf("vab %d %.20g %.20g %.20g %.20g\n", id, va_d[id], (numin - nu) * ialphaD, vb_d[id], dnu * ialphaD);
}
else{
va_d[id] = (float)(-nu * ialphaD);
vb_d[id] = (float)(ialphaD);
}
vyf_d[id] = (float)(vy_d[id]);
Sf_d[id] = (float)(S_d[id]);
S1f_d[id] = (float)(S1_d[id]);
//if(id < 100) printf("S3 %d %g %g %g\n", id, S_d[id], S1_d[id], vy_d[id]);
}
}
__global__ void L_kernelExomol(double *readBuffer_d, double *nu_d, double *S_d, double *EL_d, double *ialphaD_d, double* A_d, double *vy_d, double *n_d, const double defaultL, const double defaultn, const double gammaF, const double mass, const double T, const double Q, const double Abundance, const double Sscale, const int NL, const int kk){
int id = blockIdx.x * blockDim.x + threadIdx.x + kk;
if(id < NL){
nu_d[id] = readBuffer_d[id * 4 + 0];
double S = readBuffer_d[id * 4 + 1];
EL_d[id] = readBuffer_d[id * 4 + 2];
double A = readBuffer_d[id * 4 + 3];
//if(id < 10) printf("%d %g %g %g %g\n", id, nu_d[id], S, EL_d[id], A);
ialphaD_d[id] = def_c * sqrt( mass / (2.0 * def_kB * T));
A_d[id] = A / (4.0 * M_PI * def_c);
vy_d[id] = defaultL * gammaF;
n_d[id] = defaultn;
S_d[id] = S * Abundance * Sscale / Q;
// if(id < 100) printf("%d %g %g %g %g %g %g %g\n", id, nu_d[id], S_d[id], ialphaD_d[id], EL_d[id], Q, 0.0, vy_d[id]);
}
}
__global__ void L_kernelExomolSuper(double *readBuffer_d, double *nu_d, double *S_d, double *ialphaD_d, double *vy_d, double *n_d, const double defaultL, const double defaultn, const double gammaF, const double mass, const double T, const double Abundance, const double Sscale, const int NL, const int kk){
int id = blockIdx.x * blockDim.x + threadIdx.x + kk;
if(id < NL){
nu_d[id] = readBuffer_d[id * 2 + 0];
double S = readBuffer_d[id * 2 + 1];
//if(id < 10) printf("%d %g %g\n", id, nu_d[id], S);
ialphaD_d[id] = def_c * sqrt( mass / (2.0 * def_kB * T));
vy_d[id] = defaultL * gammaF;
n_d[id] = defaultn;
S_d[id] = S * Abundance * Sscale;
// if(id < 100) printf("%d %g %g %g %g %g %g %g\n", id, nu_d[id], S_d[id], ialphaD_d[id], EL_d[id], Q, 0.0, vy_d[id]);
}
}
__global__ void L_kernelKurucz(double *readBuffer_d, double *nu_d, double *S_d, double *EL_d, double *ialphaD_d, double* A_d, double *vy_d, double *n_d, const double defaultL, const double defaultn, const double gammaF, const double mass, const double T, const double Q, const double Abundance, const double Sscale, const int NL, const int kk){
int id = blockIdx.x * blockDim.x + threadIdx.x + kk;
if(id < NL){
nu_d[id] = readBuffer_d[id * 5 + 0];
double S = readBuffer_d[id * 5 + 1];
EL_d[id] = readBuffer_d[id * 5 + 2];
double A = readBuffer_d[id * 5 + 3];
double GammaN = readBuffer_d[id * 5 + 4];
//if(id < 10) printf("%d %g %g %g %g %g\n", id, nu_d[id], S, EL_d[id], A, GammaN);
ialphaD_d[id] = def_c * sqrt( mass / (2.0 * def_kB * T));
A_d[id] = (A + GammaN) / (4.0 * M_PI * def_c);
vy_d[id] = defaultL * gammaF;
n_d[id] = defaultn;
S_d[id] = S * Abundance * Sscale / Q;
// if(id < 100) printf("%d %g %g %g %g %g %g %g\n", id, nu_d[id], S_d[id], ialphaD_d[id], EL_d[id], Q, 0.0, vy_d[id]);
}
}
__global__ void L_kernelNIST(double *readBuffer_d, double *nu_d, double *S_d, double *EL_d, double *ialphaD_d, double* A_d, double *vy_d, double *n_d, const double defaultL, const double defaultn, const double gammaF, const double mass, const double T, const double Q, const double Abundance, const double Sscale, const int NL, const int kk){
int id = blockIdx.x * blockDim.x + threadIdx.x + kk;
if(id < NL){
nu_d[id] = readBuffer_d[id * 5 + 0];
double S = readBuffer_d[id * 5 + 1];
EL_d[id] = readBuffer_d[id * 5 + 2];
double A = readBuffer_d[id * 5 + 3];
double GammaN = readBuffer_d[id * 5 + 4];
//if(id < 10) printf("%d %g %g %g %g %g\n", id, nu_d[id], S, EL_d[id], A, GammaN);
ialphaD_d[id] = def_c * sqrt( mass / (2.0 * def_kB * T));
A_d[id] = (A + GammaN) / (4.0 * M_PI * def_c);
vy_d[id] = defaultL * gammaF;
n_d[id] = defaultn;
S_d[id] = S * Abundance * Sscale / Q;
// if(id < 100) printf("%d %g %g %g %g %g %g %g\n", id, nu_d[id], S_d[id], ialphaD_d[id], EL_d[id], Q, 0.0, vy_d[id]);
}
}
__global__ void L_kernelVALD(double *readBuffer_d, double *nu_d, double *S_d, double *EL_d, double *ialphaD_d, double* A_d, double *vy_d, double *n_d, const double defaultL, const double defaultn, const double gammaF, const double mass, const double T, const double Q, const double Abundance, const double Sscale, const int NL, const int kk){
int id = blockIdx.x * blockDim.x + threadIdx.x + kk;
if(id < NL){
nu_d[id] = readBuffer_d[id * 5 + 0];
double S = readBuffer_d[id * 5 + 1];
EL_d[id] = readBuffer_d[id * 5 + 2];
double A = readBuffer_d[id * 5 + 3];
double GammaN = readBuffer_d[id * 5 + 4];
//if(id < 10) printf("%d %g %g %g %g %g\n", id, nu_d[id], S, EL_d[id], A, GammaN);
ialphaD_d[id] = def_c * sqrt( mass / (2.0 * def_kB * T));
A_d[id] = (A + GammaN) / (4.0 * M_PI * def_c);
vy_d[id] = defaultL * gammaF;
n_d[id] = defaultn;
S_d[id] = S * Abundance * Sscale / Q;
// if(id < 100) printf("%d %g %g %g %g %g %g %g\n", id, nu_d[id], S_d[id], ialphaD_d[id], EL_d[id], Q, 0.0, vy_d[id]);
}
}
// *************************************************
//This kernel calls directly the Voigt function
// It is usefull to test the profile
//
//Author Simon Grimm
//November 2014
// *************************************************
__global__ void Voigt_line_kernel(double a, double dnu, double *K_d, double Nx, double xmax){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx;
if(id < Nx){
double aTOL = M_PI * sqrt(-1.0 / log(def_TOL * 0.5));
double x = fabs(-xmax + id * 2.0 * xmax / ((double)(Nx)));
K_d[id] = voigt_916(x, a, aTOL, id);
}
}
// *************************************************
//This kernel calls directly the Voigt function
// It is usefull to test the profile
//
//Author Simon Grimm
//November 2014
// *************************************************
__global__ void Voigt_2d_kernel(const double a, const double b, const double c, double *K_d, int Nx, int Ny, size_t pitch, double xMax, double yMax){
int idx = blockIdx.x * blockDim.x + threadIdx.x;
int idy = blockIdx.y * blockDim.y + threadIdx.y;
if(idx < Nx && idy < Ny){
//x and y arrays from 0.1 to 2000
double x = exp(-2.3 + idx * xMax / double(Nx - 1));
double y = exp(-2.3 + idy * yMax / double(Ny - 1));
double s1, s2, s3;
double ex2 = expf(-x * x);
//Compute Sigma Series
if(x != 0.0 && y != 0.0) Sigmab(x, y, s1, s2, s3, a, ex2, idx);
double xy = x * y;
double cos2xy = cosf(2.0 * xy);
double sinxy = sinf(xy);
double t1 = ex2 * erfcx(y) * cos2xy;
double t2 = sinxy * ex2 * sinxy / y;
double t3 = y * (-cos2xy * s1 + 0.5 * (s2 + s3));
t1 += c * (t2 + t3);
if(x == 0.0) t1 = erfcx(y);
if(y == 0.0) t1 = ex2;
//K_d[idy * Nx + idx] = t1 * b;
double xxyy = x * x + y * y;
double *row = (double *)(((char *)K_d)+(idy*pitch));
if(xxyy < 100.0){
row[idx] = t1 * b;
}
else if(xxyy < 1.0e6){
//Region B
double t1 = x * x * 6.0;
double t2 = xxyy + 1.5;
double t3 = M_PI * t2;
double t4 = (t3 * (2.0 * t2 + xxyy) - 2.0 * t1) / (3.0 * xxyy * (t3 * t2 - t1));
row[idx] = y * t4 / M_PI;
}
else{
//Region A
row[idx] = y / (M_PI * xxyy);
}
}
}
__global__ void Voigt_2df_kernel(const float a, const float b, const float c, float *K_d, int Nx, int Ny, size_t pitch, float xMax, float yMax){
int idx = blockIdx.x * blockDim.x + threadIdx.x;
int idy = blockIdx.y * blockDim.y + threadIdx.y;
if(idx < Nx && idy < Ny){
float x = idx * xMax / float(Nx - 1);
float y = idy * yMax / float(Ny - 1);
float s1, s2, s3;
float ex2 = expf(-x * x);
//Compute Sigma Series
if(x != 0.0f && y != 0.0f) Sigmabf(x, y, s1, s2, s3, a, ex2, idx);
float xy = x * y;
float cos2xy = cosf(2.0f * xy);
float sinxy = sinf(xy);
float t1 = ex2 * erfcxf(y) * cos2xy;
float t2 = sinxy * ex2 * sinxy / y;
float t3 = y * (-cos2xy * s1 + 0.5f * (s2 + s3));
t1 += c * (t2 + t3);
if(x == 0.0f) t1 = erfcxf(y);
if(y == 0.0f) t1 = ex2;
//K_d[idy * Nx + idx] = t1 * b;
float *row = (float *)(((char *)K_d)+(idy*pitch));
row[idx] = t1 * b;
//row[idx] = float(idy * Nx + idx);
//if(idy == 0) printf("a %d %d %g %g %g\n", idx, idy, x, y, float(idy * Nx + idx));
//printf("%g %g %g %g %g %g\n", x, y, s1, s2, s3, K_d[idy * Nx + idx]);
}
}
// *************************************************
//This kernel initializes K_d with kmin
//
//Author Simon Grimm
//September 2016
// *************************************************
__global__ void InitialK_kernel(double *K_d, const int Nx, const double kmin, const int k){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx + k;
if(id < Nx){
K_d[id] = kmin;
}
}
// *******************************************************
// This kernel adds the different streams together to K_d
// Author: Simon Grimm
// May 2019
// ********************************************************
__global__ void AddKStreams_kernel(double *K_d, double *KS_d, const int nStreams, const int NX){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx;
if(id < NX){
double K = 0.0;
for(int i = 0; i < nStreams; ++i){
K += KS_d[i * NX + id];
//if(id == 0) printf("KS %d %d %g %g\n", id, i, K, KS_d[i * NX + id]);
KS_d[i * NX + id] = 0.0;
}
K_d[id] += K;
}
}
// *************************************************
// This kernel initializes the location of nu
//
// Author Simon Grimm
// September 2016
// *************************************************
__global__ void setX_kernel(double *x_d, const double Nx, const double numin, const double dnu, const int Nxb, int useIndividualX, double *binBoundaries_d, const int k){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx + k;
if(id < Nx){
if(useIndividualX == 0){
x_d[id] = numin + id * dnu;
}
else{
int bin = id / Nxb;
double dnu = (binBoundaries_d[bin + 1] - binBoundaries_d[bin]) / ((double)(Nxb));
int start = bin * Nxb;
x_d[id] = binBoundaries_d[bin] + (id - start) * dnu;
}
//printf("x_d %d %.20g\n", id, x_d[id]);
}
}
// *************************************************
//This kernel initializes the binkeys with the bin number
//
//Nx is the total number of points,
//Nxb is the number of points ber bin
//
//Author Simon Grimm
//Septermber 2016
// *************************************************
__global__ void binKey_kernel(int *binKey_d, const int Nx, const int Nxb, double *binBoundaries_d, const int nbins, const double numax, double *x_d, int useIndividualX, const int k){
int id = blockIdx.x * blockDim.x + threadIdx.x + k;
if(id < Nx){
//int bin = id / Nxb;
//binKey_d[id] = bin;
//if(useIndividualX == 1){
double nu = x_d[id];
binKey_d[id] = 0;
for(int i = 0; i < nbins; ++i){
double nl = binBoundaries_d[i];
double nr = binBoundaries_d[i + 1];
if(nl <= nu && nu < nr){
binKey_d[id] = i;
break;
}
}
if(nu >= numax){
binKey_d[id] = nbins;
}
//}
}
}
// *************************************************
//This kernel determines the starting index of the bins
//
//Nx is the total number of points,
//nbins is the number of bins
//
//Author Simon Grimm
//September 2016
// *************************************************
__global__ void binIndex_kernel(int *binKey_d, int *binIndex_d, const int Nx, const int nbins, const int k){
int id = blockIdx.x * blockDim.x + threadIdx.x + k;
if(id < Nx - 1){
int bin = binKey_d[id];
int bin1 = binKey_d[id + 1];
if(bin1 > bin){
binIndex_d[bin1] = id + 1;
}
}
if(id == 0){
binIndex_d[0] = 0;
binIndex_d[nbins] = Nx;
binIndex_d[nbins + 1] = Nx;
}
}
//*************************************************
//This kernel copies an array a to b
//NL is the size of the array
//
//Author Simon Grimm
//January 2015
//*************************************************
__global__ void Copy_kernel(double *a_d, double *b_d, int NL, int k){
int id = blockIdx.x * blockDim.x + threadIdx.x + k;
if(id < NL){
b_d[id] = a_d[id];
}
}
//*************************************************
//This kernel sorts the array a_d acording to the key in ID_d
//and writes the result into b_d
//NL is the size of the array
//
//Author Simon Grimm
//January 2015
//*************************************************
__global__ void Sort_kernel(double *a_d, double *b_d, int *ID_d, int NL, int k){
int id = blockIdx.x * blockDim.x + threadIdx.x + k;
if(id < NL){
b_d[id] = a_d[ID_d[id]];
}
}
__global__ void print_kernel(double *nu_d, double *ialphaD_d, double *vy_d, int *ID_d, int N, int EE){
for(int id = 0; id < N; ++id){
printf("ID %d %g %g %g %d\n", EE, nu_d[id], ialphaD_d[id], vy_d[id], ID_d[id]);
}
}
// *************************************************
//This kernel computes the line shape
//
// E = 0 first order
// E = 1 third order
// E = 2 higher oder
// E = -1 first order with reduced resolution in x
// E = 10, 11 correct boundary to lower resolution
// E = 12, 13 correct boundary to cut off
//
//Author Simon Grimm
//October 2016
// *************************************************
template <int NB, const int E>
__global__ void Line2f_kernel(float *S1_d, float *vy_d, float *va_d, float *vb_d, float *vcut2_d, double *K_d, const int il, const int nstart, const int Nk, const int nl, const int useIndividualX, const int Nxb, double *binBoundaries_d, const float a, const float b, const float c, const int profile){
int idx = threadIdx.x;
int id = blockIdx.x * blockDim.x + idx;
__shared__ float S1_s[NB];
__shared__ float vy_s[NB];
__shared__ float va_s[NB];
__shared__ float vb_s[NB];
__shared__ float vcut2_s[NB];
if(idx < nl){
S1_s[idx] = S1_d[il + idx];
vy_s[idx] = vy_d[il + idx];
va_s[idx] = va_d[il + idx];
vb_s[idx] = vb_d[il + idx];
vcut2_s[idx] = vcut2_d[il + idx];
}
__syncthreads();
float x;
double dnu;
double bb;
int bstart;
if(id < Nk){
int ii = nstart + id;
double K = 0.0;
if(useIndividualX != 0){
int bin = ii / Nxb;
bb = binBoundaries_d[bin];
dnu = (binBoundaries_d[bin + 1] - bb) / ((double)(Nxb));
bstart = bin * Nxb;
//printf("a %d %d %d %g %g %d %g\n", id, ii, bin, bb, dnu, bstart, dnu);
}
for(int ill = 0; ill < nl; ++ill){
float y = vy_s[ill];
if(useIndividualX == 0){
if(E >= 0){
x = va_s[ill] + ii * vb_s[ill];
}
else{
x = va_s[ill] + ii * 10 * vb_s[ill];
}
}
else{
if(E >= 0){
x = bb * vb_s[ill] + va_s[ill] + (dnu * (ii - bstart) * vb_s[ill]);
}
else{
x = bb * vb_s[ill] + va_s[ill] + (dnu * (ii * 10 - bstart) * vb_s[ill]);
}
}
//printf("x %d %d %g %g %g\n", ill, id, x, vb_s[ill], va_s[ill]);
float t1 = x * x;
float xxyy = t1 + y * y;
if(profile == 1){
if(t1 < vcut2_s[ill]){
if(E <= 0 && xxyy >= 1.0e6f){
//1 order Gauss Hermite Quadrature
K += S1_s[ill] / xxyy;
}
if(E == 1 && xxyy >= 100.0f && xxyy < 1.0e6f){
//2nd order Gauss Hermite Quadrature
t1 *= 6.0f;
float t2 = xxyy + 1.5f;
float t3 = M_PIf * t2;
float t4 = (t3 * (2.0f * t2 + xxyy) - 2.0f * t1) / (3.0f * xxyy * (t3 * t2 - t1));
K += S1_s[ill] * t4;
}
if(E == 2 && xxyy < 100.0f){
float s1, s2, s3;
float ex2 = expf(-x * x);
//Compute Sigma Series
if(x != 0.0 && y != 0.0) Sigmabf(x, y, s1, s2, s3, a, ex2, ii);
float xy = x * y;
float cos2xy = cosf(2.0f * xy);
float sinxy = sinf(xy);
float t1 = ex2 * erfcxf(y) * cos2xy;
float t2 = sinxy * ex2 * sinxy / y;
float t3 = y * (-cos2xy * s1 + 0.5f * (s2 + s3));
t1 += c * (t2 + t3);
if(x == 0.0f) t1 = erfcxf(y);
if(y == 0.0f) t1 = ex2;
K += S1_s[ill] * t1 * b;
}
if(E == 10){ //correct interpolation
int i0 = (int)((-sqrtf((1.0e6f - y * y)) - va_s[ill]) / vb_s[ill]);
i0 = (i0 / 10) * 10;
//printf("%d %d %d %d\n", id, ii, nstart, i0);
float x0 = va_s[ill] + i0 * vb_s[ill];
float xxyy0 = x0 * x0 + y * y;
float Kc0 = S1_s[ill] / xxyy0;
float Kc = Kc0 - Kc0 / 10.0 * (ii - i0);
if(Kc >= 0.0f && Kc <= Kc0){
if(xxyy >= 1.0e6f) K += S1_s[ill] / xxyy;
K -= Kc;
}
}
if(E == 11){ //correct interpolation
int i0 = (int)((sqrtf((1.0e6f - y * y)) - va_s[ill]) / vb_s[ill]);