kernel_functions.hpp 8.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228
  1. #ifndef _SCTL_KERNEL_FUNCTIONS_HPP_
  2. #define _SCTL_KERNEL_FUNCTIONS_HPP_
  3. namespace SCTL_NAMESPACE {
  4. struct Laplace3D_FxU {
  5. template <class Real> static constexpr Real ScaleFactor() {
  6. return 1 / (4 * const_pi<Real>());
  7. }
  8. template <class Real> static void Eval(Real (&u)[1][1], const Real (&r)[3], const Real (&n)[3], void* ctx_ptr) {
  9. Real r2 = r[0]*r[0]+r[1]*r[1]+r[2]*r[2];
  10. Real rinv = (r2>1e-16 ? 1/sqrt<Real>(r2) : 0);
  11. u[0][0] = rinv;
  12. }
  13. };
  14. struct Laplace3D_DxU {
  15. template <class Real> static constexpr Real ScaleFactor() {
  16. return 1 / (4 * const_pi<Real>());
  17. }
  18. template <class Real> static void Eval(Real (&u)[1][1], const Real (&r)[3], const Real (&n)[3], void* ctx_ptr) {
  19. Real r2 = r[0]*r[0]+r[1]*r[1]+r[2]*r[2];
  20. Real rinv = (r2>1e-16 ? 1/sqrt<Real>(r2) : 0);
  21. Real rdotn = r[0]*n[0] + r[1]*n[1] + r[2]*n[2];
  22. Real rinv3 = rinv * rinv * rinv;
  23. u[0][0] = rdotn * rinv3;
  24. }
  25. };
  26. struct Laplace3D_FxdU{
  27. template <class Real> static constexpr Real ScaleFactor() {
  28. return 1 / (4 * const_pi<Real>());
  29. }
  30. template <class Real> static void Eval(Real (&u)[1][3], const Real (&r)[3], const Real (&n)[3], void* ctx_ptr) {
  31. Real r2 = r[0]*r[0]+r[1]*r[1]+r[2]*r[2];
  32. Real rinv = (r2>1e-16 ? 1/sqrt<Real>(r2) : 0);
  33. Real rinv3 = rinv * rinv * rinv;
  34. u[0][0] = -r[0] * rinv3;
  35. u[0][1] = -r[1] * rinv3;
  36. u[0][2] = -r[2] * rinv3;
  37. }
  38. };
  39. template <class uKernel> class GenericKernel {
  40. template <class Real, Integer D, Integer K0, Integer K1> static constexpr Integer get_DIM (void (*uKer)(Real (&u)[K0][K1], const Real (&r)[D], const Real (&n)[D], void* ctx_ptr)) { return D; }
  41. template <class Real, Integer D, Integer K0, Integer K1> static constexpr Integer get_KDIM0(void (*uKer)(Real (&u)[K0][K1], const Real (&r)[D], const Real (&n)[D], void* ctx_ptr)) { return K0; }
  42. template <class Real, Integer D, Integer K0, Integer K1> static constexpr Integer get_KDIM1(void (*uKer)(Real (&u)[K0][K1], const Real (&r)[D], const Real (&n)[D], void* ctx_ptr)) { return K1; }
  43. static constexpr Integer DIM = get_DIM (uKernel::template Eval<double>);
  44. static constexpr Integer KDIM0 = get_KDIM0(uKernel::template Eval<double>);
  45. static constexpr Integer KDIM1 = get_KDIM1(uKernel::template Eval<double>);
  46. public:
  47. static constexpr Integer CoordDim() {
  48. return DIM;
  49. }
  50. static constexpr Integer SrcDim() {
  51. return KDIM0;
  52. }
  53. static constexpr Integer TrgDim() {
  54. return KDIM1;
  55. }
  56. template <class Real> static constexpr Real ScaleFactor() {
  57. return uKernel::template ScaleFactor<Real>();
  58. }
  59. template <class Real, Integer DOF> void uKernelEval(Iterator<Real> u, ConstIterator<Real> r, ConstIterator<Real> n, ConstIterator<Real> f) const {
  60. Real M[KDIM0][KDIM1];
  61. uKernel::Eval(M, *(Real (*)[DIM])r, *(Real (*)[DIM])n, ctx_ptr);
  62. for (Integer i = 0; i < DOF; i++) {
  63. for (Integer k0 = 0; k0 < KDIM0; k0++) {
  64. for (Integer k1 = 0; k1 < KDIM1; k1++) {
  65. u[i*KDIM1+k1] += M[k0][k1] * f[i*KDIM0+k0];
  66. }
  67. }
  68. }
  69. }
  70. template <class Real, Integer DOF_ = 0> void Eval(Vector<Real>& U, const Vector<Real>& Xt, const Vector<Real>& Xs, const Vector<Real>& Xn, const Vector<Real>& F) const {
  71. Long Ns = Xs.Dim() / DIM;
  72. Long Nt = Xt.Dim() / DIM;
  73. Long DOF = (DOF_ ? DOF_ : F.Dim() / Ns / KDIM0);
  74. SCTL_ASSERT(Xs.Dim() == Ns * DIM);
  75. SCTL_ASSERT(Xt.Dim() == Nt * DIM);
  76. SCTL_ASSERT(F.Dim() == Ns * DOF * KDIM0);
  77. if (!DOF_) {
  78. if (DOF == 1) Eval<Real,1>(U,Xt,Xs,Xn,F);
  79. if (DOF == 2) Eval<Real,2>(U,Xt,Xs,Xn,F);
  80. if (DOF == 3) Eval<Real,3>(U,Xt,Xs,Xn,F);
  81. if (DOF == 4) Eval<Real,4>(U,Xt,Xs,Xn,F);
  82. if (DOF == 5) Eval<Real,5>(U,Xt,Xs,Xn,F);
  83. if (DOF == 6) Eval<Real,6>(U,Xt,Xs,Xn,F);
  84. if (DOF <= 6) return;
  85. }
  86. if (U.Dim() != Nt * DOF * KDIM1) {
  87. U.ReInit(Nt * DOF * KDIM1);
  88. U.SetZero();
  89. }
  90. #pragma omp parallel for schedule(static)
  91. for (Long t = 0; t < Nt; t++) {
  92. Iterator<Real> u;
  93. Vector<Real> u_buff0;
  94. StaticArray<Real,DOF_*KDIM1> u_buff1;
  95. if (!DOF_) { // Set u
  96. u_buff0.ReInit(DOF*KDIM1);
  97. u = u_buff0.begin();
  98. } else {
  99. u = (Iterator<Real>)u_buff1;
  100. }
  101. for (Integer k = 0; k < DOF*KDIM1; k++) {
  102. u[k] = 0;
  103. }
  104. StaticArray<Real,DIM> r;
  105. for (Long s = 0; s < Ns; s++) {
  106. for (Integer k = 0; k < DIM; k++) {
  107. r[k] = Xs[s*DIM+k] - Xt[t*DIM+k];
  108. }
  109. ConstIterator<Real> n = Xn.begin() + s*DIM;
  110. if (DOF_) {
  111. ConstIterator<Real> f = F.begin() + s*DOF_*KDIM0;
  112. this->template uKernelEval<Real,DOF_>(u,r,n,f);
  113. } else {
  114. ConstIterator<Real> f = F.begin() + s*DOF*KDIM0;
  115. for (Integer k = 0; k < DOF; k++) {
  116. this->template uKernelEval<Real,1>(u + k*KDIM1,r,n,f + k*KDIM0);
  117. }
  118. }
  119. }
  120. for (Integer k = 0; k < DOF*KDIM1; k++) {
  121. U[t*DOF*KDIM1 + k] += u[k];
  122. }
  123. }
  124. }
  125. template <class Real> void KernelMatrix(Matrix<Real>& M, const Vector<Real>& Xt, const Vector<Real>& Xs, const Vector<Real>& Xn) const {
  126. Long Ns = Xs.Dim() / DIM;
  127. Long Nt = Xt.Dim() / DIM;
  128. SCTL_ASSERT(Xs.Dim() == Ns * DIM);
  129. SCTL_ASSERT(Xt.Dim() == Nt * DIM);
  130. if (M.Dim(0) != Ns * KDIM0 || M.Dim(1) != Nt * KDIM1) {
  131. M.ReInit(Ns * KDIM0, Nt * KDIM1);
  132. }
  133. #pragma omp parallel for schedule(static)
  134. for (Long t = 0; t < Nt; t++) {
  135. StaticArray<Real,DIM> r;
  136. StaticArray<Real,KDIM0> f;
  137. for (Integer k = 0; k < KDIM0; k++) {
  138. f[k] = 0;
  139. }
  140. for (Long s = 0; s < Ns; s++) {
  141. for (Integer k = 0; k < DIM; k++) {
  142. r[k] = Xs[s*DIM+k] - Xt[t*DIM+k];
  143. }
  144. ConstIterator<Real> n = Xn.begin() + s*DIM;
  145. for (Integer k = 0; k < KDIM0; k++) {
  146. f[k] = 1;
  147. Iterator<Real> u = M[s*KDIM0+k] + t*KDIM1;
  148. for (Integer i = 0; i < KDIM1; i++) u[i] = 0;
  149. this->template uKernelEval<Real,1>(u,r,n,f);
  150. f[k] = 0;
  151. }
  152. }
  153. }
  154. }
  155. private:
  156. void* ctx_ptr;
  157. };
  158. template <class Real, Integer DIM> class ParticleFMM {
  159. public:
  160. template <class Kernel> static void Eval(Vector<Real>& U, const Vector<Real>& Xt, const Vector<Real>& Xs, const Vector<Real>& Xn, const Vector<Real>& F, const Kernel& kernel, const Comm& comm) {
  161. #ifdef SCTL_HAVE_PVFMM
  162. SCTL_ASSERT(false);
  163. #else
  164. constexpr Integer KDIM0 = Kernel::SrcDim();
  165. constexpr Integer KDIM1 = Kernel::TrgDim();
  166. const Integer rank = comm.Rank();
  167. const Integer np = comm.Size();
  168. const Long Ns = Xs.Dim() / DIM;
  169. const Long Nt = Xt.Dim() / DIM;
  170. SCTL_ASSERT(Xt.Dim() == Nt * DIM);
  171. const Long dof = F.Dim() / (Ns * KDIM0);
  172. SCTL_ASSERT(Ns && F.Dim() == Ns * dof * KDIM0);
  173. Vector<Real> Xs_, Xn_, F_;
  174. if (U.Dim() != Nt * dof * KDIM1) {
  175. U.ReInit(Nt * dof * KDIM1);
  176. U.SetZero();
  177. }
  178. for (Long i = 0; i < np; i++) {
  179. auto send_recv_vec = [comm,rank,np](Vector<Real>& X_, const Vector<Real>& X, Integer offset){
  180. Integer send_partner = (rank + offset) % np;
  181. Integer recv_partner = (rank + np - offset) % np;
  182. Long send_cnt = X.Dim(), recv_cnt;
  183. void* recv_req = comm.Irecv( Ptr2Itr<Long>(&recv_cnt,1), 1, recv_partner, offset);
  184. void* send_req = comm.Isend(Ptr2ConstItr<Long>(&send_cnt,1), 1, send_partner, offset);
  185. comm.Wait(recv_req);
  186. comm.Wait(send_req);
  187. X_.ReInit(recv_cnt);
  188. recv_req = comm.Irecv(X_.begin(), recv_cnt, recv_partner, offset);
  189. send_req = comm.Isend(X .begin(), send_cnt, send_partner, offset);
  190. comm.Wait(recv_req);
  191. comm.Wait(send_req);
  192. };
  193. send_recv_vec(Xs_, Xs, i);
  194. send_recv_vec(Xn_, Xn, i);
  195. send_recv_vec(F_ , F , i);
  196. kernel.Eval(U, Xt, Xs_, Xn_, F_);
  197. }
  198. #endif
  199. }
  200. private:
  201. };
  202. } // end namespace
  203. #endif //_SCTL_KERNEL_FUNCTIONS_HPP_