@@ -36,56 +36,59 @@ namespace picongpu
36
36
{
37
37
using namespace PMacc;
38
38
39
- template<class EBox>
40
- __global__ void kernelLaserE(EBox fieldE, LaserManipulator lMan)
39
+ struct kernelLaserE
41
40
{
42
- DataSpace<simDim> cellOffset;
43
-
44
- cellOffset.x() = (blockIdx.x * blockDim.x) + threadIdx.x;
45
- DataSpace<simDim> EFieldOffset;
46
- EFieldOffset.x() = cellOffset.x() + (MappingDesc::SuperCellSize::x::value * GUARD_SIZE);
47
- EFieldOffset.y() = MappingDesc::SuperCellSize::y::value*GUARD_SIZE;
48
-
49
- #if (SIMDIM==DIM3)
50
- cellOffset.z() = (blockIdx.y * blockDim.y) + threadIdx.y;
51
- EFieldOffset.z() = cellOffset.z() + (MappingDesc::SuperCellSize::z::value * GUARD_SIZE);
52
- #endif
53
-
54
- //uint32_t zOffset
55
-
56
- /** Calculate how many neighbors to the left we have
57
- * to initialize the laser in the E-Field
58
- *
59
- * Example: Yee needs one neighbor to perform dB = curlE
60
- * -> initialize in y=0 plane
61
- * A second order solver could need 2 neighbors left:
62
- * -> initialize in y=0 and y=1 plane
63
- *
64
- * Question: Why do other codes initialize the B-Field instead?
65
- * Answer: Because our fields are defined on the lower cell side
66
- * (C-Style ftw). Therefore, our curls (for example Yee)
67
- * are shifted nabla+ <-> nabla- compared to Fortran codes
68
- * (in other words: curlLeft <-> curlRight)
69
- * for E and B.
70
- * For this reason, we have to initialize E instead of B.
71
- *
72
- * Problem: that's still not our case. For example our Yee does a
73
- * dE = curlLeft(B) - therefor, we should init B, too.
74
- */
75
- //const int max_y_neighbors = Get<fieldSolver::FieldSolver::OffsetOrigin_E, 1 >::value;
76
- const int max_y_neighbors = 1;
77
-
78
- for (int totalOffsetY = 0; totalOffsetY < max_y_neighbors; ++totalOffsetY)
41
+ template<class EBox>
42
+ DINLINE void operator()(EBox fieldE, LaserManipulator lMan) const
79
43
{
80
- /** \todo Right now, the phase could be wrong ( == is cloned)
81
- * \See LaserPhysics.hpp
44
+ DataSpace<simDim> cellOffset;
45
+
46
+ cellOffset.x() = (blockIdx.x * blockDim.x) + threadIdx.x;
47
+ DataSpace<simDim> EFieldOffset;
48
+ EFieldOffset.x() = cellOffset.x() + (MappingDesc::SuperCellSize::x::value * GUARD_SIZE);
49
+ EFieldOffset.y() = MappingDesc::SuperCellSize::y::value*GUARD_SIZE;
50
+
51
+ #if (SIMDIM==DIM3)
52
+ cellOffset.z() = (blockIdx.y * blockDim.y) + threadIdx.y;
53
+ EFieldOffset.z() = cellOffset.z() + (MappingDesc::SuperCellSize::z::value * GUARD_SIZE);
54
+ #endif
55
+
56
+ //uint32_t zOffset
57
+
58
+ /** Calculate how many neighbors to the left we have
59
+ * to initialize the laser in the E-Field
82
60
*
83
- * \todo What about the B-Field in the second plane?
61
+ * Example: Yee needs one neighbor to perform dB = curlE
62
+ * -> initialize in y=0 plane
63
+ * A second order solver could need 2 neighbors left:
64
+ * -> initialize in y=0 and y=1 plane
65
+ *
66
+ * Question: Why do other codes initialize the B-Field instead?
67
+ * Answer: Because our fields are defined on the lower cell side
68
+ * (C-Style ftw). Therefore, our curls (for example Yee)
69
+ * are shifted nabla+ <-> nabla- compared to Fortran codes
70
+ * (in other words: curlLeft <-> curlRight)
71
+ * for E and B.
72
+ * For this reason, we have to initialize E instead of B.
73
+ *
74
+ * Problem: that's still not our case. For example our Yee does a
75
+ * dE = curlLeft(B) - therefor, we should init B, too.
84
76
*/
85
- cellOffset.y()=totalOffsetY;
86
- fieldE(EFieldOffset) = lMan.getManipulation(cellOffset);
77
+ //const int max_y_neighbors = Get<fieldSolver::FieldSolver::OffsetOrigin_E, 1 >::value;
78
+ const int max_y_neighbors = 1;
79
+
80
+ for (int totalOffsetY = 0; totalOffsetY < max_y_neighbors; ++totalOffsetY)
81
+ {
82
+ /** \todo Right now, the phase could be wrong ( == is cloned)
83
+ * \See LaserPhysics.hpp
84
+ *
85
+ * \todo What about the B-Field in the second plane?
86
+ */
87
+ cellOffset.y()=totalOffsetY;
88
+ fieldE(EFieldOffset) = lMan.getManipulation(cellOffset);
89
+ }
87
90
}
88
- }
91
+ };
89
92
90
93
}
91
94
0 commit comments