boundary_quadrature.hpp 69 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789
  1. #ifndef _SCTL_BOUNDARY_QUADRATURE_HPP_
  2. #define _SCTL_BOUNDARY_QUADRATURE_HPP_
  3. #include <mutex>
  4. #include <atomic>
  5. #include <tuple>
  6. namespace SCTL_NAMESPACE {
  7. template <class Real, Integer DIM, Integer ORDER> class Basis {
  8. public:
  9. using ValueType = Real;
  10. // class EvalOperator {
  11. // public:
  12. // };
  13. using EvalOpType = Matrix<ValueType>;
  14. static constexpr Long Dim() {
  15. return DIM;
  16. }
  17. static constexpr Long Size() {
  18. return pow<DIM,Long>(ORDER);
  19. }
  20. static const Matrix<ValueType>& Nodes() {
  21. static Matrix<ValueType> nodes_(DIM,Size());
  22. auto nodes_1d = [](Integer i) {
  23. return 0.5 - 0.5 * sctl::cos<ValueType>((2*i+1) * const_pi<ValueType>() / (2*ORDER));
  24. };
  25. { // Set nodes_
  26. static std::mutex mutex;
  27. static std::atomic<Integer> first_time(true);
  28. if (first_time.load(std::memory_order_relaxed)) {
  29. std::lock_guard<std::mutex> guard(mutex);
  30. if (first_time.load(std::memory_order_relaxed)) {
  31. Integer N = 1;
  32. for (Integer d = 0; d < DIM; d++) {
  33. for (Integer j = 0; j < ORDER; j++) {
  34. for (Integer i = 0; i < N; i++) {
  35. for (Integer k = 0; k < d; k++) {
  36. nodes_[k][j*N+i] = nodes_[k][i];
  37. }
  38. nodes_[d][j*N+i] = nodes_1d(j);
  39. }
  40. }
  41. N *= ORDER;
  42. }
  43. std::atomic_thread_fence(std::memory_order_seq_cst);
  44. first_time.store(false);
  45. }
  46. }
  47. }
  48. return nodes_;
  49. }
  50. static void Grad(Vector<Basis>& dX, const Vector<Basis>& X) {
  51. static Matrix<ValueType> GradOp[DIM];
  52. static std::mutex mutex;
  53. static std::atomic<Integer> first_time(true);
  54. if (first_time.load(std::memory_order_relaxed)) {
  55. std::lock_guard<std::mutex> guard(mutex);
  56. if (first_time.load(std::memory_order_relaxed)) {
  57. { // Set GradOp
  58. auto nodes = Basis<ValueType,1,ORDER>::Nodes();
  59. SCTL_ASSERT(nodes.Dim(1) == ORDER);
  60. Matrix<ValueType> M(ORDER, ORDER);
  61. for (Integer i = 0; i < ORDER; i++) { // Set M
  62. Real x = nodes[0][i];
  63. for (Integer j = 0; j < ORDER; j++) {
  64. M[j][i] = 0;
  65. for (Integer l = 0; l < ORDER; l++) {
  66. if (l != j) {
  67. Real M_ = 1;
  68. for (Integer k = 0; k < ORDER; k++) {
  69. if (k != j && k != l) M_ *= (x - nodes[0][k]);
  70. if (k != j) M_ /= (nodes[0][j] - nodes[0][k]);
  71. }
  72. M[j][i] += M_;
  73. }
  74. }
  75. }
  76. }
  77. for (Integer d = 0; d < DIM; d++) {
  78. GradOp[d].ReInit(Size(), Size());
  79. GradOp[d] = 0;
  80. Integer stride0 = sctl::pow<Integer>(ORDER, d);
  81. Integer repeat0 = sctl::pow<Integer>(ORDER, d);
  82. Integer stride1 = sctl::pow<Integer>(ORDER, d+1);
  83. Integer repeat1 = sctl::pow<Integer>(ORDER, DIM-d-1);
  84. for (Integer k1 = 0; k1 < repeat1; k1++) {
  85. for (Integer i = 0; i < ORDER; i++) {
  86. for (Integer j = 0; j < ORDER; j++) {
  87. for (Integer k0 = 0; k0 < repeat0; k0++) {
  88. GradOp[d][k1*stride1 + i*stride0 + k0][k1*stride1 + j*stride0 + k0] = M[i][j];
  89. }
  90. }
  91. }
  92. }
  93. }
  94. }
  95. std::atomic_thread_fence(std::memory_order_seq_cst);
  96. first_time.store(false);
  97. }
  98. }
  99. if (dX.Dim() != X.Dim()*DIM) dX.ReInit(X.Dim()*DIM);
  100. for (Long i = 0; i < X.Dim(); i++) {
  101. const Matrix<ValueType> Vi(1, Size(), (Iterator<ValueType>)(ConstIterator<ValueType>)X[i].NodeValues_, false);
  102. for (Integer k = 0; k < DIM; k++) {
  103. Matrix<ValueType> Vo(1, Size(), dX[i*DIM+k].NodeValues_, false);
  104. Matrix<ValueType>::GEMM(Vo, Vi, GradOp[k]);
  105. }
  106. }
  107. }
  108. static EvalOpType SetupEval(const Matrix<ValueType>& X) {
  109. Long N = X.Dim(1);
  110. SCTL_ASSERT(X.Dim(0) == DIM);
  111. Matrix<ValueType> M(Size(), N);
  112. { // Set M
  113. auto nodes = Basis<ValueType,1,ORDER>::Nodes();
  114. Integer NN = nodes.Dim(1);
  115. Matrix<ValueType> M_(NN, DIM*N);
  116. for (Long i = 0; i < DIM*N; i++) {
  117. ValueType x = X[0][i];
  118. for (Integer j = 0; j < NN; j++) {
  119. ValueType y = 1;
  120. for (Integer k = 0; k < NN; k++) {
  121. y *= (j==k ? 1 : (nodes[0][k] - x) / (nodes[0][k] - nodes[0][j]));
  122. }
  123. M_[j][i] = y;
  124. }
  125. }
  126. if (DIM == 1) {
  127. SCTL_ASSERT(M.Dim(0) == M_.Dim(0));
  128. SCTL_ASSERT(M.Dim(1) == M_.Dim(1));
  129. M = M_;
  130. } else {
  131. Integer NNN = 1;
  132. M = 1;
  133. for (Integer d = 0; d < DIM; d++) {
  134. for (Integer k = 1; k < NN; k++) {
  135. for (Integer j = 0; j < NNN; j++) {
  136. for (Long i = 0; i < N; i++) {
  137. M[k*NNN+j][i] = M[j][i] * M_[k][d*N+i];
  138. }
  139. }
  140. }
  141. { // k = 0
  142. for (Integer j = 0; j < NNN; j++) {
  143. for (Long i = 0; i < N; i++) {
  144. M[j][i] *= M_[0][d*N+i];
  145. }
  146. }
  147. }
  148. NNN *= NN;
  149. }
  150. }
  151. }
  152. return M;
  153. }
  154. static void Eval(Matrix<ValueType>& Y, const Vector<Basis>& X, const EvalOpType& M) {
  155. Long N0 = X.Dim();
  156. Long N1 = M.Dim(1);
  157. SCTL_ASSERT(M.Dim(0) == Size());
  158. if (Y.Dim(0) != N0 || Y.Dim(1) != N1) Y.ReInit(N0, N1);
  159. for (Long i = 0; i < N0; i++) {
  160. const Matrix<ValueType> X_(1,Size(),(Iterator<ValueType>)(ConstIterator<ValueType>)X[i].NodeValues_,false);
  161. Matrix<ValueType> Y_(1,N1,Y[i],false);
  162. Matrix<ValueType>::GEMM(Y_,X_,M);
  163. }
  164. }
  165. const ValueType& operator[](Long i) const {
  166. SCTL_ASSERT(i < Size());
  167. return NodeValues_[i];
  168. }
  169. ValueType& operator[](Long i) {
  170. SCTL_ASSERT(i < Size());
  171. return NodeValues_[i];
  172. }
  173. private:
  174. StaticArray<ValueType,Size()> NodeValues_;
  175. };
  176. template <Integer COORD_DIM, class Basis> class ElemList {
  177. public:
  178. using CoordBasis = Basis;
  179. using CoordType = typename CoordBasis::ValueType;
  180. static constexpr Integer CoordDim() {
  181. return COORD_DIM;
  182. }
  183. static constexpr Integer ElemDim() {
  184. return CoordBasis::Dim();
  185. }
  186. ElemList(Long Nelem = 0) {
  187. ReInit(Nelem);
  188. }
  189. void ReInit(Long Nelem = 0) {
  190. Nelem_ = Nelem;
  191. X_.ReInit(Nelem_ * COORD_DIM);
  192. }
  193. void ReInit(const Vector<CoordBasis>& X) {
  194. Nelem_ = X.Dim() / COORD_DIM;
  195. SCTL_ASSERT(X.Dim() == Nelem_ * COORD_DIM);
  196. X_ = X;
  197. }
  198. Long NElem() const {
  199. return Nelem_;
  200. }
  201. CoordBasis& operator()(Long elem, Integer dim) {
  202. SCTL_ASSERT(elem >= 0 && elem < Nelem_);
  203. SCTL_ASSERT(dim >= 0 && dim < COORD_DIM);
  204. return X_[elem*COORD_DIM+dim];
  205. }
  206. const CoordBasis& operator()(Long elem, Integer dim) const {
  207. SCTL_ASSERT(elem >= 0 && elem < Nelem_);
  208. SCTL_ASSERT(dim >= 0 && dim < COORD_DIM);
  209. return X_[elem*COORD_DIM+dim];
  210. }
  211. const Vector<CoordBasis>& ElemVector() const {
  212. return X_;
  213. }
  214. private:
  215. static_assert(CoordBasis::Dim() <= CoordDim(), "Basis dimension can not be greater than COORD_DIM.");
  216. Vector<CoordBasis> X_;
  217. Long Nelem_;
  218. mutable Vector<CoordBasis> dX_;
  219. };
  220. template <class Real> class Quadrature {
  221. static Real machine_epsilon() {
  222. Real eps=1;
  223. while(eps*(Real)0.5+(Real)1.0>1.0) eps*=0.5;
  224. return eps;
  225. }
  226. template <Integer DIM> static void DuffyQuad(Matrix<Real>& nodes, Vector<Real>& weights, const Vector<Real>& coord, Integer order, Real adapt = -1.0) {
  227. SCTL_ASSERT(coord.Dim() == DIM);
  228. static Real eps = machine_epsilon()*16;
  229. Matrix<Real> qx;
  230. Vector<Real> qw;
  231. { // Set qx, qw
  232. Vector<Real> qx0, qw0;
  233. ChebBasis<Real>::quad_rule(order, qx0, qw0);
  234. Integer N = sctl::pow<DIM,Integer>(order);
  235. qx.ReInit(DIM,N);
  236. qw.ReInit(N);
  237. qw[0] = 1;
  238. Integer N_ = 1;
  239. for (Integer d = 0; d < DIM; d++) {
  240. for (Integer j = 0; j < order; j++) {
  241. for (Integer i = 0; i < N_; i++) {
  242. for (Integer k = 0; k < d; k++) {
  243. qx[k][j*N_+i] = qx[k][i];
  244. }
  245. qx[d][j*N_+i] = qx0[j];
  246. qw[j*N_+i] = qw[i];
  247. }
  248. }
  249. for (Integer j = 0; j < order; j++) {
  250. for (Integer i = 0; i < N_; i++) {
  251. qw[j*N_+i] *= qw0[j];
  252. }
  253. }
  254. N_ *= order;
  255. }
  256. }
  257. Vector<Real> X;
  258. { // Set X
  259. StaticArray<Real,2*DIM+2> X_;
  260. X_[0] = 0;
  261. X_[1] = adapt;
  262. for (Integer i = 0; i < DIM; i++) {
  263. X_[2*i+2] = sctl::fabs<Real>(coord[i]);
  264. X_[2*i+3] = sctl::fabs<Real>(coord[i]-1);
  265. }
  266. std::sort((Iterator<Real>)X_, (Iterator<Real>)X_+2*DIM+2);
  267. X.PushBack(std::max<Real>(0, X_[2*DIM]-1));
  268. for (Integer i = 0; i < 2*DIM+2; i++) {
  269. if (X[X.Dim()-1] < X_[i]) {
  270. if (X.Dim())
  271. X.PushBack(X_[i]);
  272. }
  273. }
  274. /////////////////////////////////////////////////////////////////////////////////////////////////
  275. Vector<Real> r(1);
  276. r[0] = X[0];
  277. for (Integer i = 1; i < X.Dim(); i++) {
  278. while (r[r.Dim() - 1] > 0.0 && (order*0.5) * r[r.Dim() - 1] < X[i]) r.PushBack((order*0.5) * r[r.Dim() - 1]); // TODO
  279. r.PushBack(X[i]);
  280. }
  281. X = r;
  282. /////////////////////////////////////////////////////////////////////////////////////////////////
  283. }
  284. Vector<Real> nds, wts;
  285. for (Integer k = 0; k < X.Dim()-1; k++) {
  286. for (Integer dd = 0; dd < 2*DIM; dd++) {
  287. Integer d0 = (dd>>1);
  288. StaticArray<Real,2*DIM> range0, range1;
  289. { // Set range0, range1
  290. Integer d1 = (dd%2?1:-1);
  291. for (Integer d = 0; d < DIM; d++) {
  292. range0[d*2+0] = std::max<Real>(0,std::min<Real>(1,coord[d] - X[k] ));
  293. range0[d*2+1] = std::max<Real>(0,std::min<Real>(1,coord[d] + X[k] ));
  294. range1[d*2+0] = std::max<Real>(0,std::min<Real>(1,coord[d] - X[k+1]));
  295. range1[d*2+1] = std::max<Real>(0,std::min<Real>(1,coord[d] + X[k+1]));
  296. }
  297. range0[d0*2+0] = std::max<Real>(0,std::min<Real>(1,coord[d0] + d1*X[k+0]));
  298. range0[d0*2+1] = std::max<Real>(0,std::min<Real>(1,coord[d0] + d1*X[k+0]));
  299. range1[d0*2+0] = std::max<Real>(0,std::min<Real>(1,coord[d0] + d1*X[k+1]));
  300. range1[d0*2+1] = std::max<Real>(0,std::min<Real>(1,coord[d0] + d1*X[k+1]));
  301. }
  302. { // if volume(range0, range1) == 0 then continue
  303. Real v0 = 1, v1 = 1;
  304. for (Integer d = 0; d < DIM; d++) {
  305. if (d == d0) {
  306. v0 *= sctl::fabs<Real>(range0[d*2+0]-range1[d*2+0]);
  307. v1 *= sctl::fabs<Real>(range0[d*2+0]-range1[d*2+0]);
  308. } else {
  309. v0 *= range0[d*2+1]-range0[d*2+0];
  310. v1 *= range1[d*2+1]-range1[d*2+0];
  311. }
  312. }
  313. if (v0 < eps && v1 < eps) continue;
  314. }
  315. for (Integer i = 0; i < qx.Dim(1); i++) { // Set nds, wts
  316. Real w = qw[i];
  317. Real z = qx[d0][i];
  318. for (Integer d = 0; d < DIM; d++) {
  319. Real y = qx[d][i];
  320. nds.PushBack((range0[d*2+0]*(1-y) + range0[d*2+1]*y)*(1-z) + (range1[d*2+0]*(1-y) + range1[d*2+1]*y)*z);
  321. if (d == d0) {
  322. w *= abs(range1[d*2+0] - range0[d*2+0]);
  323. } else {
  324. w *= (range0[d*2+1] - range0[d*2+0])*(1-z) + (range1[d*2+1] - range1[d*2+0])*z;
  325. }
  326. }
  327. wts.PushBack(w);
  328. }
  329. }
  330. }
  331. nodes = Matrix<Real>(nds.Dim()/DIM,DIM,nds.begin()).Transpose();
  332. weights = wts;
  333. }
  334. template <Integer DIM> static void TensorProductGaussQuad(Matrix<Real>& nodes, Vector<Real>& weights, Integer order) {
  335. Vector<Real> coord(DIM);
  336. coord = 0;
  337. coord[0] = -10;
  338. DuffyQuad<DIM>(nodes, weights, coord, order);
  339. }
  340. template <class DensityBasis, class ElemList, class Kernel> static void SetupSingular(Matrix<Real>& M_singular, const Matrix<Real>& trg_nds, const ElemList& elem_lst, const Kernel& kernel, Integer order_singular = 10, Integer order_direct = 10) {
  341. using CoordBasis = typename ElemList::CoordBasis;
  342. using CoordEvalOpType = typename CoordBasis::EvalOpType;
  343. using DensityEvalOpType = typename DensityBasis::EvalOpType;
  344. constexpr Integer CoordDim = ElemList::CoordDim();
  345. constexpr Integer ElemDim = ElemList::ElemDim();
  346. constexpr Integer KDIM0 = Kernel::SrcDim();
  347. constexpr Integer KDIM1 = Kernel::TrgDim();
  348. const Long Nelem = elem_lst.NElem();
  349. const Integer Ntrg = trg_nds.Dim(1);
  350. SCTL_ASSERT(trg_nds.Dim(0) == ElemDim);
  351. Vector<Real> Xt;
  352. { // Set Xt
  353. auto Meval = CoordBasis::SetupEval(trg_nds);
  354. eval_basis(Xt, elem_lst.ElemVector(), ElemList::CoordDim(), trg_nds.Dim(1), Meval);
  355. }
  356. SCTL_ASSERT(Xt.Dim() == Nelem * Ntrg * CoordDim);
  357. const Vector<CoordBasis>& X = elem_lst.ElemVector();
  358. Vector<CoordBasis> dX;
  359. CoordBasis::Grad(dX, X);
  360. auto& M = M_singular;
  361. M.ReInit(Nelem * KDIM0 * DensityBasis::Size(), KDIM1 * Ntrg);
  362. #pragma omp parallel for schedule(static)
  363. for (Long i = 0; i < Ntrg; i++) { // Set M (singular)
  364. Matrix<Real> quad_nds;
  365. Vector<Real> quad_wts;
  366. { // Set quad_nds, quad_wts
  367. StaticArray<Real,ElemDim> trg_node_;
  368. for (Integer k = 0; k < ElemDim; k++) {
  369. trg_node_[k] = trg_nds[k][i];
  370. }
  371. Vector<Real> trg_node(ElemDim, trg_node_, false);
  372. DuffyQuad<ElemDim>(quad_nds, quad_wts, trg_node, order_singular);
  373. }
  374. const CoordEvalOpType CoordEvalOp = CoordBasis::SetupEval(quad_nds);
  375. Integer Nnds = quad_wts.Dim();
  376. Vector<Real> X_, dX_, Xa_, Xn_;
  377. { // Set X_, dX_
  378. eval_basis(X_, X, CoordDim, Nnds, CoordEvalOp);
  379. eval_basis(dX_, dX, CoordDim * ElemDim, Nnds, CoordEvalOp);
  380. }
  381. if (CoordDim == 3 && ElemDim == 2) { // Compute Xa_, Xn_
  382. Long N = Nelem*Nnds;
  383. Xa_.ReInit(N);
  384. Xn_.ReInit(N*CoordDim);
  385. for (Long j = 0; j < N; j++) {
  386. StaticArray<Real,CoordDim> normal;
  387. normal[0] = dX_[j*6+2]*dX_[j*6+5] - dX_[j*6+4]*dX_[j*6+3];
  388. normal[1] = dX_[j*6+4]*dX_[j*6+1] - dX_[j*6+0]*dX_[j*6+5];
  389. normal[2] = dX_[j*6+0]*dX_[j*6+3] - dX_[j*6+2]*dX_[j*6+1];
  390. Xa_[j] = sctl::sqrt<Real>(normal[0]*normal[0]+normal[1]*normal[1]+normal[2]*normal[2]);
  391. Real invXa = 1/Xa_[j];
  392. Xn_[j*3+0] = normal[0] * invXa;
  393. Xn_[j*3+1] = normal[1] * invXa;
  394. Xn_[j*3+2] = normal[2] * invXa;
  395. }
  396. }
  397. DensityEvalOpType DensityEvalOp;
  398. if (std::is_same<CoordBasis,DensityBasis>::value) {
  399. DensityEvalOp = CoordEvalOp;
  400. } else {
  401. DensityEvalOp = DensityBasis::SetupEval(quad_nds);
  402. }
  403. for (Long j = 0; j < Nelem; j++) {
  404. Matrix<Real> M__(Nnds * KDIM0, KDIM1);
  405. { // Set kernel matrix M__
  406. const Vector<Real> X0_(CoordDim, (Iterator<Real>)Xt.begin() + (j * Ntrg + i) * CoordDim, false);
  407. const Vector<Real> X__(Nnds * CoordDim, X_.begin() + j * Nnds * CoordDim, false);
  408. const Vector<Real> Xn__(Nnds * CoordDim, Xn_.begin() + j * Nnds * CoordDim, false);
  409. kernel.template KernelMatrix<Real>(M__, X0_, X__, Xn__);
  410. }
  411. for (Long k0 = 0; k0 < KDIM0; k0++) {
  412. for (Long k1 = 0; k1 < KDIM1; k1++) {
  413. for (Long l = 0; l < DensityBasis::Size(); l++) {
  414. Real M_lk = 0;
  415. for (Long n = 0; n < Nnds; n++) {
  416. Real quad_wt = Xa_[j * Nnds + n] * quad_wts[n];
  417. M_lk += DensityEvalOp[l][n] * quad_wt * M__[n*KDIM0+k0][k1];
  418. }
  419. M[(j * KDIM0 + k0) * DensityBasis::Size() + l][k1 * Ntrg + i] = M_lk;
  420. }
  421. }
  422. }
  423. }
  424. }
  425. { // Set M (subtract direct)
  426. Matrix<Real> quad_nds;
  427. Vector<Real> quad_wts;
  428. TensorProductGaussQuad<ElemDim>(quad_nds, quad_wts, order_direct);
  429. const CoordEvalOpType CoordEvalOp = CoordBasis::SetupEval(quad_nds);
  430. Integer Nnds = quad_wts.Dim();
  431. Vector<Real> X_, dX_, Xa_, Xn_;
  432. { // Set X_, dX_
  433. eval_basis(X_, X, CoordDim, Nnds, CoordEvalOp);
  434. eval_basis(dX_, dX, CoordDim * ElemDim, Nnds, CoordEvalOp);
  435. }
  436. if (CoordDim == 3 && ElemDim == 2) { // Compute Xa_, Xn_
  437. Long N = Nelem*Nnds;
  438. Xa_.ReInit(N);
  439. Xn_.ReInit(N*CoordDim);
  440. for (Long j = 0; j < N; j++) {
  441. StaticArray<Real,CoordDim> normal;
  442. normal[0] = dX_[j*6+2]*dX_[j*6+5] - dX_[j*6+4]*dX_[j*6+3];
  443. normal[1] = dX_[j*6+4]*dX_[j*6+1] - dX_[j*6+0]*dX_[j*6+5];
  444. normal[2] = dX_[j*6+0]*dX_[j*6+3] - dX_[j*6+2]*dX_[j*6+1];
  445. Xa_[j] = sctl::sqrt<Real>(normal[0]*normal[0]+normal[1]*normal[1]+normal[2]*normal[2]);
  446. Real invXa = 1/Xa_[j];
  447. Xn_[j*3+0] = normal[0] * invXa;
  448. Xn_[j*3+1] = normal[1] * invXa;
  449. Xn_[j*3+2] = normal[2] * invXa;
  450. }
  451. }
  452. DensityEvalOpType DensityEvalOp;
  453. if (std::is_same<CoordBasis,DensityBasis>::value) {
  454. DensityEvalOp = CoordEvalOp;
  455. } else {
  456. DensityEvalOp = DensityBasis::SetupEval(quad_nds);
  457. }
  458. #pragma omp parallel for schedule(static)
  459. for (Long i = 0; i < Ntrg; i++) { // Subtract direct contribution
  460. for (Long j = 0; j < Nelem; j++) {
  461. Matrix<Real> M__(Nnds * KDIM0, KDIM1);
  462. { // Set kernel matrix M__
  463. const Vector<Real> X0_(CoordDim, (Iterator<Real>)Xt.begin() + (j * Ntrg + i) * CoordDim, false);
  464. const Vector<Real> X__(Nnds * CoordDim, X_.begin() + j * Nnds * CoordDim, false);
  465. const Vector<Real> Xn__(Nnds * CoordDim, Xn_.begin() + j * Nnds * CoordDim, false);
  466. kernel.template KernelMatrix<Real>(M__, X0_, X__, Xn__);
  467. }
  468. for (Long k0 = 0; k0 < KDIM0; k0++) {
  469. for (Long k1 = 0; k1 < KDIM1; k1++) {
  470. for (Long l = 0; l < DensityBasis::Size(); l++) {
  471. Real M_lk = 0;
  472. for (Long n = 0; n < Nnds; n++) {
  473. Real quad_wt = Xa_[j * Nnds + n] * quad_wts[n];
  474. M_lk += DensityEvalOp[l][n] * quad_wt * M__[n*KDIM0+k0][k1];
  475. }
  476. M[(j * KDIM0 + k0) * DensityBasis::Size() + l][k1 * Ntrg + i] -= M_lk;
  477. }
  478. }
  479. }
  480. }
  481. }
  482. }
  483. }
  484. template <class DensityBasis> static void EvalSingular(Matrix<Real>& U, const Vector<DensityBasis>& density, const Matrix<Real>& M, Integer KDIM0_, Integer KDIM1_) {
  485. if (M.Dim(0) == 0 || M.Dim(1) == 0) {
  486. U.ReInit(0,0);
  487. return;
  488. }
  489. const Long Ntrg = M.Dim(1) / KDIM1_;
  490. SCTL_ASSERT(M.Dim(1) == KDIM1_ * Ntrg);
  491. const Long Nelem = M.Dim(0) / (KDIM0_ * DensityBasis::Size());
  492. SCTL_ASSERT(M.Dim(0) == Nelem * KDIM0_ * DensityBasis::Size());
  493. const Integer dof = density.Dim() / (Nelem * KDIM0_);
  494. SCTL_ASSERT(density.Dim() == Nelem * dof * KDIM0_);
  495. if (U.Dim(0) != Nelem * dof * KDIM1_ || U.Dim(1) != Ntrg) {
  496. U.ReInit(Nelem * dof * KDIM1_, Ntrg);
  497. U = 0;
  498. }
  499. for (Long j = 0; j < Nelem; j++) {
  500. const Matrix<Real> M_(KDIM0_ * DensityBasis::Size(), KDIM1_ * Ntrg, (Iterator<Real>)M[j * KDIM0_ * DensityBasis::Size()], false);
  501. Matrix<Real> U_(dof, KDIM1_ * Ntrg, U[j*dof*KDIM1_], false);
  502. Matrix<Real> F_(dof, KDIM0_ * DensityBasis::Size());
  503. for (Long i = 0; i < dof; i++) {
  504. for (Long k = 0; k < KDIM0_; k++) {
  505. for (Long l = 0; l < DensityBasis::Size(); l++) {
  506. F_[i][k * DensityBasis::Size() + l] = density[(j * dof + i) * KDIM0_ + k][l];
  507. }
  508. }
  509. }
  510. Matrix<Real>::GEMM(U_, F_, M_);
  511. }
  512. }
  513. template <Integer DIM> struct PointData {
  514. bool operator<(const PointData& p) const {
  515. return mid < p.mid;
  516. }
  517. Long rank;
  518. Long surf_rank;
  519. Morton<DIM> mid;
  520. StaticArray<Real,DIM> coord;
  521. Real radius2;
  522. };
  523. template <class T1, class T2> struct Pair {
  524. Pair() {}
  525. Pair(T1 x, T2 y) : first(x), second(y) {}
  526. bool operator<(const Pair& p) const {
  527. return (first < p.first) || (((first == p.first) && (second < p.second)));
  528. }
  529. T1 first;
  530. T2 second;
  531. };
  532. template <class ElemList> static void BuildNbrList(Vector<Pair<Long,Long>>& pair_lst, const Vector<Real>& Xt, const Vector<Long>& trg_surf, const ElemList& elem_lst, Real distance_factor, Real period_length, const Comm& comm) {
  533. using CoordBasis = typename ElemList::CoordBasis;
  534. constexpr Integer CoordDim = ElemList::CoordDim();
  535. constexpr Integer ElemDim = ElemList::ElemDim();
  536. using PtData = PointData<CoordDim>;
  537. const Integer rank = comm.Rank();
  538. Real R0 = 0;
  539. StaticArray<Real,CoordDim> X0;
  540. { // Find bounding box
  541. Long N = Xt.Dim() / CoordDim;
  542. SCTL_ASSERT(Xt.Dim() == N * CoordDim);
  543. SCTL_ASSERT(N);
  544. StaticArray<Real,CoordDim*2> Xloc;
  545. StaticArray<Real,CoordDim*2> Xglb;
  546. for (Integer k = 0; k < CoordDim; k++) {
  547. Xloc[0*CoordDim+k] = Xt[k];
  548. Xloc[1*CoordDim+k] = Xt[k];
  549. }
  550. for (Long i = 0; i < N; i++) {
  551. for (Integer k = 0; k < CoordDim; k++) {
  552. Xloc[0*CoordDim+k] = std::min<Real>(Xloc[0*CoordDim+k], Xt[i*CoordDim+k]);
  553. Xloc[1*CoordDim+k] = std::max<Real>(Xloc[1*CoordDim+k], Xt[i*CoordDim+k]);
  554. }
  555. }
  556. comm.Allreduce((ConstIterator<Real>)Xloc+0*CoordDim, (Iterator<Real>)Xglb+0*CoordDim, CoordDim, Comm::CommOp::MIN);
  557. comm.Allreduce((ConstIterator<Real>)Xloc+1*CoordDim, (Iterator<Real>)Xglb+1*CoordDim, CoordDim, Comm::CommOp::MAX);
  558. for (Integer k = 0; k < CoordDim; k++) {
  559. R0 = std::max(R0, Xglb[1*CoordDim+k]-Xglb[0*CoordDim+k]);
  560. }
  561. R0 = R0 * 2.0;
  562. for (Integer k = 0; k < CoordDim; k++) {
  563. X0[k] = Xglb[k] - R0*0.25;
  564. }
  565. }
  566. if (period_length > 0) {
  567. R0 = period_length;
  568. }
  569. Vector<PtData> PtSrc, PtTrg;
  570. Integer order_upsample = (Integer)(const_pi<Real>() / distance_factor + 0.5);
  571. { // Set PtSrc
  572. const Vector<CoordBasis>& X_elem_lst = elem_lst.ElemVector();
  573. Vector<CoordBasis> dX_elem_lst;
  574. CoordBasis::Grad(dX_elem_lst, X_elem_lst);
  575. Matrix<Real> nds;
  576. Vector<Real> wts;
  577. TensorProductGaussQuad<ElemDim>(nds, wts, order_upsample);
  578. const Long Nnds = nds.Dim(1);
  579. Vector<Real> X, dX;
  580. const auto CoordEvalOp = CoordBasis::SetupEval(nds);
  581. eval_basis(X, X_elem_lst, CoordDim, Nnds, CoordEvalOp);
  582. eval_basis(dX, dX_elem_lst, CoordDim * ElemDim, Nnds, CoordEvalOp);
  583. const Long N = X.Dim() / CoordDim;
  584. const Long Nelem = elem_lst.NElem();
  585. SCTL_ASSERT(X.Dim() == N * CoordDim);
  586. SCTL_ASSERT(N == Nelem * Nnds);
  587. Long rank_offset, surf_rank_offset;
  588. { // Set rank_offset, surf_rank_offset
  589. comm.Scan(Ptr2ConstItr<Long>(&N,1), Ptr2Itr<Long>(&rank_offset,1), 1, Comm::CommOp::SUM);
  590. comm.Scan(Ptr2ConstItr<Long>(&Nelem,1), Ptr2Itr<Long>(&surf_rank_offset,1), 1, Comm::CommOp::SUM);
  591. surf_rank_offset -= Nelem;
  592. rank_offset -= N;
  593. }
  594. PtSrc.ReInit(N);
  595. const Real R0inv = 1.0 / R0;
  596. for (Long i = 0; i < N; i++) { // Set coord
  597. for (Integer k = 0; k < CoordDim; k++) {
  598. PtSrc[i].coord[k] = (X[i*CoordDim+k] - X0[k]) * R0inv;
  599. }
  600. }
  601. if (period_length > 0) { // Wrap-around coord
  602. for (Long i = 0; i < N; i++) {
  603. auto& x = PtSrc[i].coord;
  604. for (Integer k = 0; k < CoordDim; k++) {
  605. x[k] -= (Long)(x[k]);
  606. }
  607. }
  608. }
  609. for (Long i = 0; i < N; i++) { // Set radius2, mid, rank
  610. Integer depth = 0;
  611. { // Set radius2, depth
  612. Real radius2 = 0;
  613. for (Integer k0 = 0; k0 < ElemDim; k0++) {
  614. Real R2 = 0;
  615. for (Integer k1 = 0; k1 < CoordDim; k1++) {
  616. Real dX_ = dX[(i*CoordDim+k1)*ElemDim+k0];
  617. R2 += dX_*dX_;
  618. }
  619. radius2 = std::max(radius2, R2);
  620. }
  621. radius2 *= R0inv*R0inv * distance_factor*distance_factor;
  622. PtSrc[i].radius2 = radius2;
  623. Long Rinv = (Long)(1.0/radius2);
  624. while (Rinv > 0) {
  625. Rinv = (Rinv>>2);
  626. depth++;
  627. }
  628. }
  629. PtSrc[i].mid = Morton<CoordDim>((Iterator<Real>)PtSrc[i].coord, std::min(Morton<CoordDim>::MaxDepth(),depth));
  630. PtSrc[i].rank = rank_offset + i;
  631. }
  632. for (Long i = 0 ; i < Nelem; i++) { // Set surf_rank
  633. for (Long j = 0; j < Nnds; j++) {
  634. PtSrc[i*Nnds+j].surf_rank = surf_rank_offset + i;
  635. }
  636. }
  637. Vector<PtData> PtSrcSorted;
  638. comm.HyperQuickSort(PtSrc, PtSrcSorted);
  639. PtSrc.Swap(PtSrcSorted);
  640. }
  641. { // Set PtTrg
  642. const Long N = Xt.Dim() / CoordDim;
  643. SCTL_ASSERT(Xt.Dim() == N * CoordDim);
  644. Long rank_offset;
  645. { // Set rank_offset
  646. comm.Scan(Ptr2ConstItr<Long>(&N,1), Ptr2Itr<Long>(&rank_offset,1), 1, Comm::CommOp::SUM);
  647. rank_offset -= N;
  648. }
  649. PtTrg.ReInit(N);
  650. const Real R0inv = 1.0 / R0;
  651. for (Long i = 0; i < N; i++) { // Set coord
  652. for (Integer k = 0; k < CoordDim; k++) {
  653. PtTrg[i].coord[k] = (Xt[i*CoordDim+k] - X0[k]) * R0inv;
  654. }
  655. }
  656. if (period_length > 0) { // Wrap-around coord
  657. for (Long i = 0; i < N; i++) {
  658. auto& x = PtTrg[i].coord;
  659. for (Integer k = 0; k < CoordDim; k++) {
  660. x[k] -= (Long)(x[k]);
  661. }
  662. }
  663. }
  664. for (Long i = 0; i < N; i++) { // Set radius2, mid, rank
  665. PtTrg[i].radius2 = 0;
  666. PtTrg[i].mid = Morton<CoordDim>((Iterator<Real>)PtTrg[i].coord);
  667. PtTrg[i].rank = rank_offset + i;
  668. }
  669. if (trg_surf.Dim()) { // Set surf_rank
  670. SCTL_ASSERT(trg_surf.Dim() == N);
  671. for (Long i = 0; i < N; i++) {
  672. PtTrg[i].surf_rank = trg_surf[i];
  673. }
  674. } else {
  675. for (Long i = 0; i < N; i++) {
  676. PtTrg[i].surf_rank = -1;
  677. }
  678. }
  679. Vector<PtData> PtTrgSorted;
  680. comm.HyperQuickSort(PtTrg, PtTrgSorted);
  681. PtTrg.Swap(PtTrgSorted);
  682. }
  683. Tree<CoordDim> tree(comm);
  684. { // Init tree
  685. Vector<Real> Xall(PtSrc.Dim()+PtTrg.Dim());
  686. { // Set Xall
  687. Xall.ReInit((PtSrc.Dim()+PtTrg.Dim())*CoordDim);
  688. Long Nsrc = PtSrc.Dim();
  689. Long Ntrg = PtTrg.Dim();
  690. for (Long i = 0; i < Nsrc; i++) {
  691. for (Integer k = 0; k < CoordDim; k++) {
  692. Xall[i*CoordDim+k] = PtSrc[i].coord[k];
  693. }
  694. }
  695. for (Long i = 0; i < Ntrg; i++) {
  696. for (Integer k = 0; k < CoordDim; k++) {
  697. Xall[(Nsrc+i)*CoordDim+k] = PtTrg[i].coord[k];
  698. }
  699. }
  700. }
  701. tree.UpdateRefinement(Xall, 1000, true, period_length>0);
  702. }
  703. { // Repartition PtSrc, PtTrg
  704. PtData splitter;
  705. splitter.mid = tree.GetPartitionMID()[rank];
  706. comm.PartitionS(PtSrc, splitter);
  707. comm.PartitionS(PtTrg, splitter);
  708. }
  709. { // Add tree data PtSrc
  710. const auto& node_mid = tree.GetNodeMID();
  711. const Long N = node_mid.Dim();
  712. SCTL_ASSERT(N);
  713. Vector<Long> dsp(N), cnt(N);
  714. for (Long i = 0; i < N; i++) {
  715. PtData m0;
  716. m0.mid = node_mid[i];
  717. dsp[i] = std::lower_bound(PtSrc.begin(), PtSrc.end(), m0) - PtSrc.begin();
  718. }
  719. for (Long i = 0; i < N-1; i++) {
  720. cnt[i] = dsp[i+1] - dsp[i];
  721. }
  722. cnt[N-1] = PtSrc.Dim() - dsp[N-1];
  723. tree.AddData("PtSrc", PtSrc, cnt);
  724. }
  725. tree.template Broadcast<PtData>("PtSrc");
  726. { // Build pair_lst
  727. Vector<Long> cnt;
  728. Vector<PtData> PtSrc;
  729. tree.GetData(PtSrc, cnt, "PtSrc");
  730. const auto& node_mid = tree.GetNodeMID();
  731. const auto& node_attr = tree.GetNodeAttr();
  732. Vector<Morton<CoordDim>> nbr_mid_tmp;
  733. for (Long i = 0; i < node_mid.Dim(); i++) {
  734. if (node_attr[i].Leaf && !node_attr[i].Ghost) {
  735. Vector<Morton<CoordDim>> child_mid;
  736. node_mid[i].Children(child_mid);
  737. for (const auto& trg_mid : child_mid) {
  738. Integer d0 = trg_mid.Depth();
  739. Vector<PtData> Src, Trg;
  740. { // Set Trg
  741. PtData m0, m1;
  742. m0.mid = trg_mid;
  743. m1.mid = trg_mid.Next();
  744. Long a = std::lower_bound(PtTrg.begin(), PtTrg.end(), m0) - PtTrg.begin();
  745. Long b = std::lower_bound(PtTrg.begin(), PtTrg.end(), m1) - PtTrg.begin();
  746. Trg.ReInit(b-a, PtTrg.begin()+a, false);
  747. if (!Trg.Dim()) continue;
  748. }
  749. Vector<std::set<Long>> near_elem(Trg.Dim());
  750. for (Integer d = 0; d <= d0; d++) {
  751. trg_mid.NbrList(nbr_mid_tmp, d, period_length>0);
  752. for (const auto& src_mid : nbr_mid_tmp) { // Set Src
  753. PtData m0, m1;
  754. m0.mid = src_mid;
  755. m1.mid = (d==d0 ? src_mid.Next() : src_mid.Ancestor(d+1));
  756. Long a = std::lower_bound(PtSrc.begin(), PtSrc.end(), m0) - PtSrc.begin();
  757. Long b = std::lower_bound(PtSrc.begin(), PtSrc.end(), m1) - PtSrc.begin();
  758. Src.ReInit(b-a, PtSrc.begin()+a, false);
  759. if (!Src.Dim()) continue;
  760. for (Long t = 0; t < Trg.Dim(); t++) { // set near_elem[t] <-- {s : dist(s,t) < radius(s)}
  761. for (Long s = 0; s < Src.Dim(); s++) {
  762. if (Trg[t].surf_rank != Src[s].surf_rank) {
  763. Real R2 = 0;
  764. for (Integer k = 0; k < CoordDim; k++) {
  765. Real dx = (Src[s].coord[k] - Trg[t].coord[k]);
  766. R2 += dx * dx;
  767. }
  768. if (R2 < Src[s].radius2) {
  769. near_elem[t].insert(Src[s].surf_rank);
  770. }
  771. }
  772. }
  773. }
  774. }
  775. }
  776. for (Long t = 0; t < Trg.Dim(); t++) { // Set pair_lst
  777. for (Long elem_idx : near_elem[t]) {
  778. pair_lst.PushBack(Pair<Long,Long>(elem_idx,Trg[t].rank));
  779. }
  780. }
  781. }
  782. }
  783. }
  784. }
  785. { // Sort and repartition pair_lst
  786. Vector<Pair<Long,Long>> pair_lst_sorted;
  787. comm.HyperQuickSort(pair_lst, pair_lst_sorted);
  788. Long surf_rank_offset;
  789. const Long Nelem = elem_lst.NElem();
  790. comm.Scan(Ptr2ConstItr<Long>(&Nelem,1), Ptr2Itr<Long>(&surf_rank_offset,1), 1, Comm::CommOp::SUM);
  791. surf_rank_offset -= Nelem;
  792. comm.PartitionS(pair_lst_sorted, Pair<Long,Long>(surf_rank_offset,0));
  793. pair_lst.Swap(pair_lst_sorted);
  794. }
  795. }
  796. template <class ElemList> static void BuildNbrListDeprecated(Vector<Pair<Long,Long>>& pair_lst, const Vector<Real>& Xt, const ElemList& elem_lst, const Matrix<Real>& surf_nds, Real distance_factor) {
  797. using CoordBasis = typename ElemList::CoordBasis;
  798. constexpr Integer CoordDim = ElemList::CoordDim();
  799. constexpr Integer ElemDim = ElemList::ElemDim();
  800. const Long Nelem = elem_lst.NElem();
  801. const Long Ntrg = Xt.Dim() / CoordDim;
  802. SCTL_ASSERT(Xt.Dim() == Ntrg * CoordDim);
  803. Long Nnds, Nsurf_nds;
  804. Vector<Real> X_surf, X, dX;
  805. Integer order_upsample = (Integer)(const_pi<Real>() / distance_factor + 0.5);
  806. { // Set X, dX
  807. const Vector<CoordBasis>& X_elem_lst = elem_lst.ElemVector();
  808. Vector<CoordBasis> dX_elem_lst;
  809. CoordBasis::Grad(dX_elem_lst, X_elem_lst);
  810. Matrix<Real> nds_upsample;
  811. Vector<Real> wts_upsample;
  812. TensorProductGaussQuad<ElemDim>(nds_upsample, wts_upsample, order_upsample);
  813. Nnds = nds_upsample.Dim(1);
  814. const auto CoordEvalOp = CoordBasis::SetupEval(nds_upsample);
  815. eval_basis(X, X_elem_lst, CoordDim, nds_upsample.Dim(1), CoordEvalOp);
  816. eval_basis(dX, dX_elem_lst, CoordDim * ElemDim, nds_upsample.Dim(1), CoordEvalOp);
  817. Nsurf_nds = surf_nds.Dim(1);
  818. const auto CoordEvalOp_surf = CoordBasis::SetupEval(surf_nds);
  819. eval_basis(X_surf, X_elem_lst, CoordDim, Nsurf_nds, CoordEvalOp_surf);
  820. }
  821. Real d2 = distance_factor * distance_factor;
  822. for (Long i = 0; i < Nelem; i++) {
  823. std::set<Long> near_pts;
  824. std::set<Long> self_pts;
  825. for (Long j = 0; j < Nnds; j++) {
  826. Real R2_max = 0;
  827. StaticArray<Real, CoordDim> X0;
  828. for (Integer k = 0; k < CoordDim; k++) {
  829. X0[k] = X[(i*Nnds+j)*CoordDim+k];
  830. }
  831. for (Integer k0 = 0; k0 < ElemDim; k0++) {
  832. Real R2 = 0;
  833. for (Integer k1 = 0; k1 < CoordDim; k1++) {
  834. Real dX_ = dX[((i*Nnds+j)*CoordDim+k1)*ElemDim+k0];
  835. R2 += dX_*dX_;
  836. }
  837. R2_max = std::max(R2_max, R2*d2);
  838. }
  839. for (Long k = 0; k < Ntrg; k++) {
  840. Real R2 = 0;
  841. for (Integer l = 0; l < CoordDim; l++) {
  842. Real dX = Xt[k*CoordDim+l]- X0[l];
  843. R2 += dX * dX;
  844. }
  845. if (R2 < R2_max) near_pts.insert(k);
  846. }
  847. }
  848. for (Long j = 0; j < Nsurf_nds; j++) {
  849. StaticArray<Real, CoordDim> X0;
  850. for (Integer k = 0; k < CoordDim; k++) {
  851. X0[k] = X_surf[(i*Nsurf_nds+j)*CoordDim+k];
  852. }
  853. for (Long k = 0; k < Ntrg; k++) {
  854. Real R2 = 0;
  855. for (Integer l = 0; l < CoordDim; l++) {
  856. Real dX = Xt[k*CoordDim+l]- X0[l];
  857. R2 += dX * dX;
  858. }
  859. if (R2 == 0) self_pts.insert(k);
  860. }
  861. }
  862. for (Long trg_idx : self_pts) {
  863. near_pts.erase(trg_idx);
  864. }
  865. for (Long trg_idx : near_pts) {
  866. pair_lst.PushBack(Pair<Long,Long>(i,trg_idx));
  867. }
  868. }
  869. }
  870. template <class DensityBasis, class ElemList, class Kernel> static void SetupNearSingular(Matrix<Real>& M_near_singular, Vector<Pair<Long,Long>>& pair_lst, const Vector<Real>& Xt_, const Vector<Long>& trg_surf, const ElemList& elem_lst, const Kernel& kernel, Integer order_singular, Integer order_direct, Real period_length, const Comm& comm) {
  871. static_assert(std::is_same<Real,typename DensityBasis::ValueType>::value);
  872. static_assert(std::is_same<Real,typename ElemList::CoordType>::value);
  873. static_assert(DensityBasis::Dim() == ElemList::ElemDim());
  874. using CoordBasis = typename ElemList::CoordBasis;
  875. using CoordEvalOpType = typename CoordBasis::EvalOpType;
  876. using DensityEvalOpType = typename DensityBasis::EvalOpType;
  877. constexpr Integer CoordDim = ElemList::CoordDim();
  878. constexpr Integer ElemDim = ElemList::ElemDim();
  879. constexpr Integer KDIM0 = Kernel::SrcDim();
  880. constexpr Integer KDIM1 = Kernel::TrgDim();
  881. const Long Nelem = elem_lst.NElem();
  882. BuildNbrList(pair_lst, Xt_, trg_surf, elem_lst, 2.5/order_direct, period_length, comm);
  883. const Long Ninterac = pair_lst.Dim();
  884. Vector<Real> Xt;
  885. { // Set Xt
  886. Integer rank = comm.Rank();
  887. Integer np = comm.Size();
  888. Vector<Long> splitter_ranks;
  889. { // Set splitter_ranks
  890. Vector<Long> cnt(np);
  891. const Long N = Xt_.Dim() / CoordDim;
  892. comm.Allgather(Ptr2ConstItr<Long>(&N,1), 1, cnt.begin(), 1);
  893. scan(splitter_ranks, cnt);
  894. }
  895. Vector<Long> scatter_index, recv_index, recv_cnt(np), recv_dsp(np);
  896. { // Set scatter_index, recv_index, recv_cnt, recv_dsp
  897. { // Set scatter_index, recv_index
  898. Vector<Pair<Long,Long>> scatter_pair(pair_lst.Dim());
  899. for (Long i = 0; i < pair_lst.Dim(); i++) {
  900. scatter_pair[i] = Pair<Long,Long>(pair_lst[i].second,i);
  901. }
  902. omp_par::merge_sort(scatter_pair.begin(), scatter_pair.end());
  903. recv_index.ReInit(scatter_pair.Dim());
  904. scatter_index.ReInit(scatter_pair.Dim());
  905. for (Long i = 0; i < scatter_index.Dim(); i++) {
  906. recv_index[i] = scatter_pair[i].first;
  907. scatter_index[i] = scatter_pair[i].second;
  908. }
  909. }
  910. for (Integer i = 0; i < np; i++) {
  911. recv_dsp[i] = std::lower_bound(recv_index.begin(), recv_index.end(), splitter_ranks[i]) - recv_index.begin();
  912. }
  913. for (Integer i = 0; i < np-1; i++) {
  914. recv_cnt[i] = recv_dsp[i+1] - recv_dsp[i];
  915. }
  916. recv_cnt[np-1] = recv_index.Dim() - recv_dsp[np-1];
  917. }
  918. Vector<Long> send_index, send_cnt(np), send_dsp(np);
  919. { // Set send_index, send_cnt, send_dsp
  920. comm.Alltoall(recv_cnt.begin(), 1, send_cnt.begin(), 1);
  921. scan(send_dsp, send_cnt);
  922. send_index.ReInit(send_cnt[np-1] + send_dsp[np-1]);
  923. comm.Alltoallv(recv_index.begin(), recv_cnt.begin(), recv_dsp.begin(), send_index.begin(), send_cnt.begin(), send_dsp.begin());
  924. }
  925. Vector<Real> Xt_send(send_index.Dim() * CoordDim);
  926. for (Long i = 0; i < send_index.Dim(); i++) { // Set Xt_send
  927. Long idx = send_index[i] - splitter_ranks[rank];
  928. for (Integer k = 0; k < CoordDim; k++) {
  929. Xt_send[i*CoordDim+k] = Xt_[idx*CoordDim+k];
  930. }
  931. }
  932. Vector<Real> Xt_recv(recv_index.Dim() * CoordDim);
  933. { // Set Xt_recv
  934. for (Long i = 0; i < np; i++) {
  935. send_cnt[i] *= CoordDim;
  936. send_dsp[i] *= CoordDim;
  937. recv_cnt[i] *= CoordDim;
  938. recv_dsp[i] *= CoordDim;
  939. }
  940. comm.Alltoallv(Xt_send.begin(), send_cnt.begin(), send_dsp.begin(), Xt_recv.begin(), recv_cnt.begin(), recv_dsp.begin());
  941. }
  942. Xt.ReInit(scatter_index.Dim() * CoordDim);
  943. for (Long i = 0; i < scatter_index.Dim(); i++) { // Set Xt
  944. Long idx = scatter_index[i];
  945. for (Integer k = 0; k < CoordDim; k++) {
  946. Xt[idx*CoordDim+k] = Xt_recv[i*CoordDim+k];
  947. }
  948. }
  949. }
  950. const Vector<CoordBasis>& X = elem_lst.ElemVector();
  951. Vector<CoordBasis> dX;
  952. CoordBasis::Grad(dX, X);
  953. Long elem_rank_offset;
  954. { // Set elem_rank_offset
  955. comm.Scan(Ptr2ConstItr<Long>(&Nelem,1), Ptr2Itr<Long>(&elem_rank_offset,1), 1, Comm::CommOp::SUM);
  956. elem_rank_offset -= Nelem;
  957. }
  958. auto& M = M_near_singular;
  959. M.ReInit(Ninterac * KDIM0 * DensityBasis::Size(), KDIM1);
  960. #pragma omp parallel for schedule(static)
  961. for (Long j = 0; j < Ninterac; j++) { // Set M (near-singular)
  962. const Long src_idx = pair_lst[j].first - elem_rank_offset;
  963. Real adapt = -1.0;
  964. Tensor<Real,true,ElemDim,1> u0;
  965. { // Set u0 (project target point to the surface patch in parameter space)
  966. ConstIterator<Real> Xt_ = Xt.begin() + j * CoordDim;
  967. const auto& nodes = CoordBasis::Nodes();
  968. Long min_idx = -1;
  969. Real min_R2 = 1e10;
  970. for (Long i = 0; i < CoordBasis::Size(); i++) {
  971. Real R2 = 0;
  972. for (Integer k = 0; k < CoordDim; k++) {
  973. Real dX = X[src_idx * CoordDim + k][i] - Xt_[k];
  974. R2 += dX * dX;
  975. }
  976. if (R2 < min_R2) {
  977. min_R2 = R2;
  978. min_idx = i;
  979. }
  980. }
  981. SCTL_ASSERT(min_idx >= 0);
  982. for (Integer k = 0; k < ElemDim; k++) {
  983. u0(k,0) = nodes[k][min_idx];
  984. }
  985. for (Integer i = 0; i < 2; i++) { // iterate
  986. Matrix<Real> X_, dX_;
  987. for (Integer k = 0; k < ElemDim; k++) {
  988. u0(k,0) = std::min(1.0, u0(k,0));
  989. u0(k,0) = std::max(0.0, u0(k,0));
  990. }
  991. const auto eval_op = CoordBasis::SetupEval(Matrix<Real>(ElemDim,1,u0.begin(),false));
  992. CoordBasis::Eval(X_, Vector<CoordBasis>(CoordDim,(Iterator<CoordBasis>)X.begin()+src_idx*CoordDim,false),eval_op);
  993. CoordBasis::Eval(dX_, Vector<CoordBasis>(CoordDim*ElemDim,dX.begin()+src_idx*CoordDim*ElemDim,false),eval_op);
  994. const Tensor<Real,false,CoordDim,1> x0((Iterator<Real>)Xt_);
  995. const Tensor<Real,false,CoordDim,1> x(X_.begin());
  996. const Tensor<Real,false,CoordDim,ElemDim> x_u(dX_.begin());
  997. auto inv = [](const Tensor<Real,true,2,2>& M) {
  998. Tensor<Real,true,2,2> Minv;
  999. Real det_inv = 1.0 / (M(0,0)*M(1,1) - M(1,0)*M(0,1));
  1000. Minv(0,0) = M(1,1) * det_inv;
  1001. Minv(0,1) =-M(0,1) * det_inv;
  1002. Minv(1,0) =-M(1,0) * det_inv;
  1003. Minv(1,1) = M(0,0) * det_inv;
  1004. return Minv;
  1005. };
  1006. auto du = inv(x_u.RotateRight()*x_u) * x_u.RotateRight()*(x0-x);
  1007. u0 = u0 + du;
  1008. auto x_u_squared = x_u.RotateRight() * x_u;
  1009. adapt = sctl::sqrt<Real>( ((x0-x).RotateRight()*(x0-x))(0,0) / std::max<Real>(x_u_squared(0,0),x_u_squared(1,1)) );
  1010. }
  1011. }
  1012. Matrix<Real> quad_nds;
  1013. Vector<Real> quad_wts;
  1014. DuffyQuad<ElemDim>(quad_nds, quad_wts, Vector<Real>(ElemDim,u0.begin(),false), order_singular, adapt);
  1015. const CoordEvalOpType CoordEvalOp = CoordBasis::SetupEval(quad_nds);
  1016. Integer Nnds = quad_wts.Dim();
  1017. Vector<Real> X_, dX_, Xa_, Xn_;
  1018. { // Set X_, dX_
  1019. const Vector<CoordBasis> X__(CoordDim, (Iterator<CoordBasis>)X.begin() + src_idx * CoordDim, false);
  1020. const Vector<CoordBasis> dX__(CoordDim * ElemDim, (Iterator<CoordBasis>)dX.begin() + src_idx * CoordDim * ElemDim, false);
  1021. eval_basis(X_, X__, CoordDim, Nnds, CoordEvalOp);
  1022. eval_basis(dX_, dX__, CoordDim * ElemDim, Nnds, CoordEvalOp);
  1023. }
  1024. if (CoordDim == 3 && ElemDim == 2) { // Compute Xa_, Xn_
  1025. Xa_.ReInit(Nnds);
  1026. Xn_.ReInit(Nnds*CoordDim);
  1027. for (Long j = 0; j < Nnds; j++) {
  1028. StaticArray<Real,CoordDim> normal;
  1029. normal[0] = dX_[j*6+2]*dX_[j*6+5] - dX_[j*6+4]*dX_[j*6+3];
  1030. normal[1] = dX_[j*6+4]*dX_[j*6+1] - dX_[j*6+0]*dX_[j*6+5];
  1031. normal[2] = dX_[j*6+0]*dX_[j*6+3] - dX_[j*6+2]*dX_[j*6+1];
  1032. Xa_[j] = sctl::sqrt<Real>(normal[0]*normal[0]+normal[1]*normal[1]+normal[2]*normal[2]);
  1033. Real invXa = 1/Xa_[j];
  1034. Xn_[j*3+0] = normal[0] * invXa;
  1035. Xn_[j*3+1] = normal[1] * invXa;
  1036. Xn_[j*3+2] = normal[2] * invXa;
  1037. }
  1038. }
  1039. DensityEvalOpType DensityEvalOp;
  1040. if (std::is_same<CoordBasis,DensityBasis>::value) {
  1041. DensityEvalOp = CoordEvalOp;
  1042. } else {
  1043. DensityEvalOp = DensityBasis::SetupEval(quad_nds);
  1044. }
  1045. Matrix<Real> M__(Nnds * KDIM0, KDIM1);
  1046. { // Set kernel matrix M__
  1047. const Vector<Real> X0_(CoordDim, (Iterator<Real>)Xt.begin() + j * CoordDim, false);
  1048. kernel.template KernelMatrix<Real>(M__, X0_, X_, Xn_);
  1049. }
  1050. for (Long k0 = 0; k0 < KDIM0; k0++) {
  1051. for (Long k1 = 0; k1 < KDIM1; k1++) {
  1052. for (Long l = 0; l < DensityBasis::Size(); l++) {
  1053. Real M_lk = 0;
  1054. for (Long n = 0; n < Nnds; n++) {
  1055. Real quad_wt = Xa_[n] * quad_wts[n];
  1056. M_lk += DensityEvalOp[l][n] * quad_wt * M__[n*KDIM0+k0][k1];
  1057. }
  1058. M[(j * KDIM0 + k0) * DensityBasis::Size() + l][k1] = M_lk;
  1059. }
  1060. }
  1061. }
  1062. }
  1063. { // Set M (subtract direct)
  1064. Matrix<Real> quad_nds;
  1065. Vector<Real> quad_wts;
  1066. TensorProductGaussQuad<ElemDim>(quad_nds, quad_wts, order_direct);
  1067. const CoordEvalOpType CoordEvalOp = CoordBasis::SetupEval(quad_nds);
  1068. Integer Nnds = quad_wts.Dim();
  1069. Vector<Real> X_, dX_, Xa_, Xn_;
  1070. { // Set X_, dX_
  1071. eval_basis(X_, X, CoordDim, Nnds, CoordEvalOp);
  1072. eval_basis(dX_, dX, CoordDim * ElemDim, Nnds, CoordEvalOp);
  1073. }
  1074. if (CoordDim == 3 && ElemDim == 2) { // Compute Xa_, Xn_
  1075. Long N = Nelem*Nnds;
  1076. Xa_.ReInit(N);
  1077. Xn_.ReInit(N*CoordDim);
  1078. for (Long j = 0; j < N; j++) {
  1079. StaticArray<Real,CoordDim> normal;
  1080. normal[0] = dX_[j*6+2]*dX_[j*6+5] - dX_[j*6+4]*dX_[j*6+3];
  1081. normal[1] = dX_[j*6+4]*dX_[j*6+1] - dX_[j*6+0]*dX_[j*6+5];
  1082. normal[2] = dX_[j*6+0]*dX_[j*6+3] - dX_[j*6+2]*dX_[j*6+1];
  1083. Xa_[j] = sctl::sqrt<Real>(normal[0]*normal[0]+normal[1]*normal[1]+normal[2]*normal[2]);
  1084. Real invXa = 1/Xa_[j];
  1085. Xn_[j*3+0] = normal[0] * invXa;
  1086. Xn_[j*3+1] = normal[1] * invXa;
  1087. Xn_[j*3+2] = normal[2] * invXa;
  1088. }
  1089. }
  1090. DensityEvalOpType DensityEvalOp;
  1091. if (std::is_same<CoordBasis,DensityBasis>::value) {
  1092. DensityEvalOp = CoordEvalOp;
  1093. } else {
  1094. DensityEvalOp = DensityBasis::SetupEval(quad_nds);
  1095. }
  1096. #pragma omp parallel for schedule(static)
  1097. for (Long j = 0; j < Ninterac; j++) { // Subtract direct contribution
  1098. const Long src_idx = pair_lst[j].first - elem_rank_offset;
  1099. Matrix<Real> M__(Nnds * KDIM0, KDIM1);
  1100. { // Set kernel matrix M__
  1101. const Vector<Real> X0_(CoordDim, (Iterator<Real>)Xt.begin() + j * CoordDim, false);
  1102. Vector<Real> X__(Nnds * CoordDim, X_.begin() + src_idx * Nnds * CoordDim, false);
  1103. Vector<Real> Xn__(Nnds * CoordDim, Xn_.begin() + src_idx * Nnds * CoordDim, false);
  1104. kernel.template KernelMatrix<Real>(M__, X0_, X__, Xn__);
  1105. }
  1106. for (Long k0 = 0; k0 < KDIM0; k0++) {
  1107. for (Long k1 = 0; k1 < KDIM1; k1++) {
  1108. for (Long l = 0; l < DensityBasis::Size(); l++) {
  1109. Real M_lk = 0;
  1110. for (Long n = 0; n < Nnds; n++) {
  1111. Real quad_wt = Xa_[src_idx * Nnds + n] * quad_wts[n];
  1112. M_lk += DensityEvalOp[l][n] * quad_wt * M__[n*KDIM0+k0][k1];
  1113. }
  1114. M[(j * KDIM0 + k0) * DensityBasis::Size() + l][k1] -= M_lk;
  1115. }
  1116. }
  1117. }
  1118. }
  1119. }
  1120. }
  1121. template <class DensityBasis> static void EvalNearSingular(Vector<Real>& U, const Vector<DensityBasis>& density, const Matrix<Real>& M, const Vector<Pair<Long,Long>>& pair_lst, Long Nelem_, Long Ntrg_, Integer KDIM0_, Integer KDIM1_, const Comm& comm) {
  1122. const Long Ninterac = pair_lst.Dim();
  1123. const Integer dof = density.Dim() / Nelem_ / KDIM0_;
  1124. SCTL_ASSERT(density.Dim() == Nelem_ * dof * KDIM0_);
  1125. Long elem_rank_offset;
  1126. { // Set elem_rank_offset
  1127. comm.Scan(Ptr2ConstItr<Long>(&Nelem_,1), Ptr2Itr<Long>(&elem_rank_offset,1), 1, Comm::CommOp::SUM);
  1128. elem_rank_offset -= Nelem_;
  1129. }
  1130. Vector<Real> U_loc(Ninterac*dof*KDIM1_);
  1131. for (Long j = 0; j < Ninterac; j++) {
  1132. const Long src_idx = pair_lst[j].first - elem_rank_offset;
  1133. const Matrix<Real> M_(KDIM0_ * DensityBasis::Size(), KDIM1_, (Iterator<Real>)M[j * KDIM0_ * DensityBasis::Size()], false);
  1134. Matrix<Real> U_(dof, KDIM1_, U_loc.begin() + j*dof*KDIM1_, false);
  1135. Matrix<Real> F_(dof, KDIM0_ * DensityBasis::Size());
  1136. for (Long i = 0; i < dof; i++) {
  1137. for (Long k = 0; k < KDIM0_; k++) {
  1138. for (Long l = 0; l < DensityBasis::Size(); l++) {
  1139. F_[i][k * DensityBasis::Size() + l] = density[(src_idx * dof + i) * KDIM0_ + k][l];
  1140. }
  1141. }
  1142. }
  1143. Matrix<Real>::GEMM(U_, F_, M_);
  1144. }
  1145. if (U.Dim() != Ntrg_ * dof * KDIM1_) {
  1146. U.ReInit(Ntrg_ * dof * KDIM1_);
  1147. U = 0;
  1148. }
  1149. { // Set U
  1150. Integer rank = comm.Rank();
  1151. Integer np = comm.Size();
  1152. Vector<Long> splitter_ranks;
  1153. { // Set splitter_ranks
  1154. Vector<Long> cnt(np);
  1155. comm.Allgather(Ptr2ConstItr<Long>(&Ntrg_,1), 1, cnt.begin(), 1);
  1156. scan(splitter_ranks, cnt);
  1157. }
  1158. Vector<Long> scatter_index, send_index, send_cnt(np), send_dsp(np);
  1159. { // Set scatter_index, send_index, send_cnt, send_dsp
  1160. { // Set scatter_index, send_index
  1161. Vector<Pair<Long,Long>> scatter_pair(pair_lst.Dim());
  1162. for (Long i = 0; i < pair_lst.Dim(); i++) {
  1163. scatter_pair[i] = Pair<Long,Long>(pair_lst[i].second,i);
  1164. }
  1165. omp_par::merge_sort(scatter_pair.begin(), scatter_pair.end());
  1166. send_index.ReInit(scatter_pair.Dim());
  1167. scatter_index.ReInit(scatter_pair.Dim());
  1168. for (Long i = 0; i < scatter_index.Dim(); i++) {
  1169. send_index[i] = scatter_pair[i].first;
  1170. scatter_index[i] = scatter_pair[i].second;
  1171. }
  1172. }
  1173. for (Integer i = 0; i < np; i++) {
  1174. send_dsp[i] = std::lower_bound(send_index.begin(), send_index.end(), splitter_ranks[i]) - send_index.begin();
  1175. }
  1176. for (Integer i = 0; i < np-1; i++) {
  1177. send_cnt[i] = send_dsp[i+1] - send_dsp[i];
  1178. }
  1179. send_cnt[np-1] = send_index.Dim() - send_dsp[np-1];
  1180. }
  1181. Vector<Long> recv_index, recv_cnt(np), recv_dsp(np);
  1182. { // Set recv_index, recv_cnt, recv_dsp
  1183. comm.Alltoall(send_cnt.begin(), 1, recv_cnt.begin(), 1);
  1184. scan(recv_dsp, recv_cnt);
  1185. recv_index.ReInit(recv_cnt[np-1] + recv_dsp[np-1]);
  1186. comm.Alltoallv(send_index.begin(), send_cnt.begin(), send_dsp.begin(), recv_index.begin(), recv_cnt.begin(), recv_dsp.begin());
  1187. }
  1188. Vector<Real> U_send(scatter_index.Dim() * dof * KDIM1_);
  1189. for (Long i = 0; i < scatter_index.Dim(); i++) {
  1190. Long idx = scatter_index[i]*dof*KDIM1_;
  1191. for (Long k = 0; k < dof * KDIM1_; k++) {
  1192. U_send[i*dof*KDIM1_ + k] = U_loc[idx + k];
  1193. }
  1194. }
  1195. Vector<Real> U_recv(recv_index.Dim() * dof * KDIM1_);
  1196. { // Set U_recv
  1197. for (Long i = 0; i < np; i++) {
  1198. send_cnt[i] *= dof * KDIM1_;
  1199. send_dsp[i] *= dof * KDIM1_;
  1200. recv_cnt[i] *= dof * KDIM1_;
  1201. recv_dsp[i] *= dof * KDIM1_;
  1202. }
  1203. comm.Alltoallv(U_send.begin(), send_cnt.begin(), send_dsp.begin(), U_recv.begin(), recv_cnt.begin(), recv_dsp.begin());
  1204. }
  1205. for (Long i = 0; i < recv_index.Dim(); i++) { // Set U
  1206. Long idx = (recv_index[i] - splitter_ranks[rank]) * dof * KDIM1_;
  1207. for (Integer k = 0; k < dof * KDIM1_; k++) {
  1208. U[idx + k] += U_recv[i*dof*KDIM1_ + k];
  1209. }
  1210. }
  1211. }
  1212. }
  1213. template <class ElemList, class DensityBasis, class Kernel> static void Direct(Vector<Real>& U, const Vector<Real>& Xt, const ElemList& elem_lst, const Vector<DensityBasis>& density, const Kernel& kernel, Integer order_direct, const Comm& comm) {
  1214. using CoordBasis = typename ElemList::CoordBasis;
  1215. using CoordEvalOpType = typename CoordBasis::EvalOpType;
  1216. using DensityEvalOpType = typename DensityBasis::EvalOpType;
  1217. constexpr Integer CoordDim = ElemList::CoordDim();
  1218. constexpr Integer ElemDim = ElemList::ElemDim();
  1219. constexpr Integer KDIM0 = Kernel::SrcDim();
  1220. constexpr Integer KDIM1 = Kernel::TrgDim();
  1221. const Long Nelem = elem_lst.NElem();
  1222. const Integer dof = density.Dim() / Nelem / KDIM0;
  1223. SCTL_ASSERT(density.Dim() == Nelem * dof * KDIM0);
  1224. Matrix<Real> quad_nds;
  1225. Vector<Real> quad_wts;
  1226. TensorProductGaussQuad<ElemDim>(quad_nds, quad_wts, order_direct);
  1227. const CoordEvalOpType CoordEvalOp = CoordBasis::SetupEval(quad_nds);
  1228. Integer Nnds = quad_wts.Dim();
  1229. const Vector<CoordBasis>& X = elem_lst.ElemVector();
  1230. Vector<CoordBasis> dX;
  1231. CoordBasis::Grad(dX, X);
  1232. Vector<Real> X_, dX_, Xa_, Xn_;
  1233. eval_basis(X_, X, CoordDim, Nnds, CoordEvalOp);
  1234. eval_basis(dX_, dX, CoordDim*ElemDim, Nnds, CoordEvalOp);
  1235. if (CoordDim == 3 && ElemDim == 2) { // Compute Xa_, Xn_
  1236. Long N = Nelem*Nnds;
  1237. Xa_.ReInit(N);
  1238. Xn_.ReInit(N*CoordDim);
  1239. for (Long j = 0; j < N; j++) {
  1240. StaticArray<Real,CoordDim> normal;
  1241. normal[0] = dX_[j*6+2]*dX_[j*6+5] - dX_[j*6+4]*dX_[j*6+3];
  1242. normal[1] = dX_[j*6+4]*dX_[j*6+1] - dX_[j*6+0]*dX_[j*6+5];
  1243. normal[2] = dX_[j*6+0]*dX_[j*6+3] - dX_[j*6+2]*dX_[j*6+1];
  1244. Xa_[j] = sctl::sqrt<Real>(normal[0]*normal[0]+normal[1]*normal[1]+normal[2]*normal[2]);
  1245. Real invXa = 1/Xa_[j];
  1246. Xn_[j*3+0] = normal[0] * invXa;
  1247. Xn_[j*3+1] = normal[1] * invXa;
  1248. Xn_[j*3+2] = normal[2] * invXa;
  1249. }
  1250. }
  1251. Vector<Real> Fa_;
  1252. { // Set Fa_
  1253. Vector<Real> F_;
  1254. if (std::is_same<CoordBasis,DensityBasis>::value) {
  1255. eval_basis(F_, density, dof * KDIM0, Nnds, CoordEvalOp);
  1256. } else {
  1257. const DensityEvalOpType EvalOp = DensityBasis::SetupEval(quad_nds);
  1258. eval_basis(F_, density, dof * KDIM0, Nnds, EvalOp);
  1259. }
  1260. Fa_.ReInit(F_.Dim());
  1261. const Integer DensityDOF = dof * KDIM0;
  1262. SCTL_ASSERT(F_.Dim() == Nelem * Nnds * DensityDOF);
  1263. for (Long j = 0; j < Nelem; j++) {
  1264. for (Integer k = 0; k < Nnds; k++) {
  1265. Long idx = j * Nnds + k;
  1266. Real quad_wt = Xa_[idx] * quad_wts[k];
  1267. for (Integer l = 0; l < DensityDOF; l++) {
  1268. Fa_[idx * DensityDOF + l] = F_[idx * DensityDOF + l] * quad_wt;
  1269. }
  1270. }
  1271. }
  1272. }
  1273. { // Evaluate potential
  1274. const Long Ntrg = Xt.Dim() / CoordDim;
  1275. SCTL_ASSERT(Xt.Dim() == Ntrg * CoordDim);
  1276. if (U.Dim() != Ntrg * dof * KDIM1) {
  1277. U.ReInit(Ntrg * dof * KDIM1);
  1278. U = 0;
  1279. }
  1280. ParticleFMM<Real,CoordDim>::Eval(U, Xt, X_, Xn_, Fa_, kernel, comm);
  1281. }
  1282. }
  1283. public:
  1284. template <class DensityBasis, class ElemList, class Kernel> void Setup(const ElemList& elem_lst, const Vector<Real>& Xt, const Kernel& kernel, Integer order_singular, Integer order_direct, Real period_length, const Comm& comm) {
  1285. order_direct_ = order_direct;
  1286. period_length_ = period_length;
  1287. comm_ = comm;
  1288. Profile::Tic("Setup", &comm_);
  1289. static_assert(std::is_same<Real,typename DensityBasis::ValueType>::value);
  1290. static_assert(std::is_same<Real,typename ElemList::CoordType>::value);
  1291. static_assert(DensityBasis::Dim() == ElemList::ElemDim());
  1292. Xt_ = Xt;
  1293. M_singular.ReInit(0,0);
  1294. Profile::Tic("SetupNearSingular", &comm_);
  1295. SetupNearSingular<DensityBasis>(M_near_singular, pair_lst, Xt_, Vector<Long>(), elem_lst, kernel, order_singular, order_direct_, period_length_, comm_);
  1296. Profile::Toc();
  1297. Profile::Toc();
  1298. }
  1299. template <class DensityBasis, class PotentialBasis, class ElemList, class Kernel> void Setup(const ElemList& elem_lst, const Kernel& kernel, Integer order_singular, Integer order_direct, Real period_length, const Comm& comm) {
  1300. order_direct_ = order_direct;
  1301. period_length_ = period_length;
  1302. comm_ = comm;
  1303. Profile::Tic("Setup", &comm_);
  1304. static_assert(std::is_same<Real,typename PotentialBasis::ValueType>::value);
  1305. static_assert(std::is_same<Real,typename DensityBasis::ValueType>::value);
  1306. static_assert(std::is_same<Real,typename ElemList::CoordType>::value);
  1307. static_assert(PotentialBasis::Dim() == ElemList::ElemDim());
  1308. static_assert(DensityBasis::Dim() == ElemList::ElemDim());
  1309. Vector<Long> trg_surf;
  1310. { // Set Xt_
  1311. using CoordBasis = typename ElemList::CoordBasis;
  1312. Matrix<Real> trg_nds = PotentialBasis::Nodes();
  1313. auto Meval = CoordBasis::SetupEval(trg_nds);
  1314. eval_basis(Xt_, elem_lst.ElemVector(), ElemList::CoordDim(), trg_nds.Dim(1), Meval);
  1315. { // Set trg_surf
  1316. const Long Nelem = elem_lst.NElem();
  1317. const Long Nnds = trg_nds.Dim(1);
  1318. Long elem_offset;
  1319. { // Set elem_offset
  1320. comm.Scan(Ptr2ConstItr<Long>(&Nelem,1), Ptr2Itr<Long>(&elem_offset,1), 1, Comm::CommOp::SUM);
  1321. elem_offset -= Nelem;
  1322. }
  1323. trg_surf.ReInit(elem_lst.NElem() * trg_nds.Dim(1));
  1324. for (Long i = 0; i < Nelem; i++) {
  1325. for (Long j = 0; j < Nnds; j++) {
  1326. trg_surf[i*Nnds+j] = elem_offset + i;
  1327. }
  1328. }
  1329. }
  1330. }
  1331. Profile::Tic("SetupSingular", &comm_);
  1332. SetupSingular<DensityBasis>(M_singular, PotentialBasis::Nodes(), elem_lst, kernel, order_singular, order_direct_);
  1333. Profile::Toc();
  1334. Profile::Tic("SetupNearSingular", &comm_);
  1335. SetupNearSingular<DensityBasis>(M_near_singular, pair_lst, Xt_, trg_surf, elem_lst, kernel, order_singular, order_direct_, period_length_, comm_);
  1336. Profile::Toc();
  1337. Profile::Toc();
  1338. }
  1339. template <class DensityBasis, class PotentialBasis, class ElemList, class Kernel> void Eval(Vector<PotentialBasis>& U, const ElemList& elements, const Vector<DensityBasis>& F, const Kernel& kernel) {
  1340. Profile::Tic("Eval", &comm_);
  1341. Matrix<Real> U_singular;
  1342. Vector<Real> U_direct, U_near_sing;
  1343. Profile::Tic("EvalDirect", &comm_);
  1344. Direct(U_direct, Xt_, elements, F, kernel, order_direct_, comm_);
  1345. Profile::Toc();
  1346. Profile::Tic("EvalSingular", &comm_);
  1347. EvalSingular(U_singular, F, M_singular, kernel.SrcDim(), kernel.TrgDim());
  1348. Profile::Toc();
  1349. Profile::Tic("EvalNearSingular", &comm_);
  1350. EvalNearSingular(U_near_sing, F, M_near_singular, pair_lst, elements.NElem(), Xt_.Dim() / ElemList::CoordDim(), kernel.SrcDim(), kernel.TrgDim(), comm_);
  1351. SCTL_ASSERT(U_near_sing.Dim() == U_direct.Dim());
  1352. Profile::Toc();
  1353. if (U.Dim() != elements.NElem() * kernel.TrgDim()) {
  1354. U.ReInit(elements.NElem() * kernel.TrgDim());
  1355. }
  1356. for (int i = 0; i < elements.NElem(); i++) {
  1357. for (int j = 0; j < PotentialBasis::Size(); j++) {
  1358. for (int k = 0; k < kernel.TrgDim(); k++) {
  1359. Real& U_ = U[i*kernel.TrgDim()+k][j];
  1360. U_ = 0;
  1361. U_ += U_direct [(i*PotentialBasis::Size()+j)*kernel.TrgDim()+k];
  1362. U_ += U_near_sing[(i*PotentialBasis::Size()+j)*kernel.TrgDim()+k];
  1363. U_ += U_singular[i*kernel.TrgDim()+k][j];
  1364. U_ *= kernel.template ScaleFactor<Real>();
  1365. }
  1366. }
  1367. }
  1368. Profile::Toc();
  1369. }
  1370. template <class DensityBasis, class ElemList, class Kernel> void Eval(Vector<Real>& U, const ElemList& elements, const Vector<DensityBasis>& F, const Kernel& kernel) {
  1371. Profile::Tic("Eval", &comm_);
  1372. Matrix<Real> U_singular;
  1373. Vector<Real> U_direct, U_near_sing;
  1374. Profile::Tic("EvalDirect", &comm_);
  1375. Direct(U_direct, Xt_, elements, F, kernel, order_direct_, comm_);
  1376. Profile::Toc();
  1377. Profile::Tic("EvalSingular", &comm_);
  1378. EvalSingular(U_singular, F, M_singular, kernel.SrcDim(), kernel.TrgDim());
  1379. Profile::Toc();
  1380. Profile::Tic("EvalNearSingular", &comm_);
  1381. EvalNearSingular(U_near_sing, F, M_near_singular, pair_lst, elements.NElem(), Xt_.Dim() / ElemList::CoordDim(), kernel.SrcDim(), kernel.TrgDim(), comm_);
  1382. SCTL_ASSERT(U_near_sing.Dim() == U_direct.Dim());
  1383. Profile::Toc();
  1384. if (U.Dim() != U_direct.Dim()) {
  1385. U.ReInit(U_direct.Dim());
  1386. }
  1387. for (int i = 0; i < U.Dim(); i++) {
  1388. U[i] = (U_direct[i] + U_near_sing[i]) * kernel.template ScaleFactor<Real>();
  1389. }
  1390. if (U_singular.Dim(1)) {
  1391. for (int i = 0; i < elements.NElem(); i++) {
  1392. for (int j = 0; j < U_singular.Dim(1); j++) {
  1393. for (int k = 0; k < kernel.TrgDim(); k++) {
  1394. Real& U_ = U[(i*U_singular.Dim(1)+j)*kernel.TrgDim()+k];
  1395. U_ += U_singular[i*kernel.TrgDim()+k][j] * kernel.template ScaleFactor<Real>();
  1396. }
  1397. }
  1398. }
  1399. }
  1400. Profile::Toc();
  1401. }
  1402. template <Integer ORDER = 5> static void test(Integer order_singular = 10, Integer order_direct = 5, const Comm& comm = Comm::World()) {
  1403. constexpr Integer COORD_DIM = 3;
  1404. constexpr Integer ELEM_DIM = COORD_DIM-1;
  1405. using ElemList = ElemList<COORD_DIM, Basis<Real, ELEM_DIM, ORDER>>;
  1406. using DensityBasis = Basis<Real, ELEM_DIM, ORDER>;
  1407. using PotentialBasis = Basis<Real, ELEM_DIM, ORDER>;
  1408. int np = comm.Size();
  1409. int rank = comm.Rank();
  1410. auto build_torus = [rank,np](ElemList& elements, long Nt, long Np, Real Rmajor, Real Rminor){
  1411. auto nodes = ElemList::CoordBasis::Nodes();
  1412. auto torus = [](Real theta, Real phi, Real Rmajor, Real Rminor) {
  1413. Real R = Rmajor + Rminor * cos<Real>(phi);
  1414. Real X = R * cos<Real>(theta);
  1415. Real Y = R * sin<Real>(theta);
  1416. Real Z = Rminor * sin<Real>(phi);
  1417. return std::make_tuple(X,Y,Z);
  1418. };
  1419. long start = Nt*Np*(rank+0)/np;
  1420. long end = Nt*Np*(rank+1)/np;
  1421. elements.ReInit(end - start);
  1422. for (long ii = start; ii < end; ii++) {
  1423. long i = ii / Np;
  1424. long j = ii % Np;
  1425. for (int k = 0; k < nodes.Dim(1); k++) {
  1426. Real X, Y, Z;
  1427. Real theta = 2 * const_pi<Real>() * (i + nodes[0][k]) / Nt;
  1428. Real phi = 2 * const_pi<Real>() * (j + nodes[1][k]) / Np;
  1429. std::tie(X,Y,Z) = torus(theta, phi, Rmajor, Rminor);
  1430. elements(ii-start,0)[k] = X;
  1431. elements(ii-start,1)[k] = Y;
  1432. elements(ii-start,2)[k] = Z;
  1433. }
  1434. }
  1435. };
  1436. ElemList elements_src, elements_trg;
  1437. build_torus(elements_src, 28, 16, 2, 1.0);
  1438. build_torus(elements_trg, 29, 17, 2, 0.99);
  1439. Vector<Real> Xt;
  1440. Vector<PotentialBasis> U_onsurf, U_offsurf;
  1441. Vector<DensityBasis> density_sl, density_dl;
  1442. { // Set Xt, elements_src, elements_trg, density_sl, density_dl, U
  1443. Real X0[COORD_DIM] = {3,2,1};
  1444. std::function<void(Real*,Real*,Real*)> potential = [X0](Real* U, Real* X, Real* Xn) {
  1445. Real dX[COORD_DIM] = {X[0]-X0[0],X[1]-X0[1],X[2]-X0[2]};
  1446. Real Rinv = 1/sqrt(dX[0]*dX[0]+dX[1]*dX[1]+dX[2]*dX[2]);
  1447. U[0] = Rinv;
  1448. };
  1449. std::function<void(Real*,Real*,Real*)> potential_normal_derivative = [X0](Real* U, Real* X, Real* Xn) {
  1450. Real dX[COORD_DIM] = {X[0]-X0[0],X[1]-X0[1],X[2]-X0[2]};
  1451. Real Rinv = 1/sqrt(dX[0]*dX[0]+dX[1]*dX[1]+dX[2]*dX[2]);
  1452. Real RdotN = dX[0]*Xn[0]+dX[1]*Xn[1]+dX[2]*Xn[2];
  1453. U[0] = -RdotN * Rinv*Rinv*Rinv;
  1454. };
  1455. DiscretizeSurfaceFn<COORD_DIM,1>(density_sl, elements_src, potential_normal_derivative);
  1456. DiscretizeSurfaceFn<COORD_DIM,1>(density_dl, elements_src, potential);
  1457. DiscretizeSurfaceFn<COORD_DIM,1>(U_onsurf , elements_src, potential);
  1458. DiscretizeSurfaceFn<COORD_DIM,1>(U_offsurf , elements_trg, potential);
  1459. for (long i = 0; i < elements_trg.NElem(); i++) { // Set Xt
  1460. for (long j = 0; j < PotentialBasis::Size(); j++) {
  1461. for (int k = 0; k < COORD_DIM; k++) {
  1462. Xt.PushBack(elements_trg(i,k)[j]);
  1463. }
  1464. }
  1465. }
  1466. }
  1467. GenericKernel<Laplace3D_DxU> Laplace_DxU;
  1468. GenericKernel<Laplace3D_FxU> Laplace_FxU;
  1469. Profile::Enable(true);
  1470. if (1) { // Greeen's identity test (Laplace, on-surface)
  1471. Profile::Tic("OnSurface", &comm);
  1472. Quadrature<Real> quadrature_DxU, quadrature_FxU;
  1473. quadrature_FxU.Setup<DensityBasis, PotentialBasis>(elements_src, Laplace_FxU, order_singular, order_direct, -1.0, comm);
  1474. quadrature_DxU.Setup<DensityBasis, PotentialBasis>(elements_src, Laplace_DxU, order_singular, order_direct, -1.0, comm);
  1475. Vector<PotentialBasis> U_sl, U_dl;
  1476. quadrature_FxU.Eval(U_sl, elements_src, density_sl, Laplace_FxU);
  1477. quadrature_DxU.Eval(U_dl, elements_src, density_dl, Laplace_DxU);
  1478. Profile::Toc();
  1479. Real max_err = 0;
  1480. Vector<PotentialBasis> err(U_onsurf.Dim());
  1481. for (long i = 0; i < U_sl.Dim(); i++) {
  1482. for (long j = 0; j < PotentialBasis::Size(); j++) {
  1483. err[i][j] = 0.5*U_onsurf[i][j] - (U_sl[i][j] + U_dl[i][j]);
  1484. max_err = std::max<Real>(max_err, fabs(err[i][j]));
  1485. }
  1486. }
  1487. { // Print error
  1488. Real glb_err;
  1489. comm.Allreduce(Ptr2ConstItr<Real>(&max_err,1), Ptr2Itr<Real>(&glb_err,1), 1, Comm::CommOp::MAX);
  1490. if (!comm.Rank()) std::cout<<"Error = "<<glb_err<<'\n';
  1491. }
  1492. { // Write VTK output
  1493. VTUData vtu;
  1494. vtu.AddElems(elements_src, err, ORDER);
  1495. vtu.WriteVTK("err", comm);
  1496. }
  1497. { // Write VTK output
  1498. VTUData vtu;
  1499. vtu.AddElems(elements_src, U_onsurf, ORDER);
  1500. vtu.WriteVTK("U", comm);
  1501. }
  1502. }
  1503. if (1) { // Greeen's identity test (Laplace, off-surface)
  1504. Profile::Tic("OffSurface", &comm);
  1505. Quadrature<Real> quadrature_DxU, quadrature_FxU;
  1506. quadrature_FxU.Setup<DensityBasis>(elements_src, Xt, Laplace_FxU, order_singular, order_direct, -1.0, comm);
  1507. quadrature_DxU.Setup<DensityBasis>(elements_src, Xt, Laplace_DxU, order_singular, order_direct, -1.0, comm);
  1508. Vector<Real> U_sl, U_dl;
  1509. quadrature_FxU.Eval(U_sl, elements_src, density_sl, Laplace_FxU);
  1510. quadrature_DxU.Eval(U_dl, elements_src, density_dl, Laplace_DxU);
  1511. Profile::Toc();
  1512. Real max_err = 0;
  1513. Vector<PotentialBasis> err(elements_trg.NElem());
  1514. for (long i = 0; i < elements_trg.NElem(); i++) {
  1515. for (long j = 0; j < PotentialBasis::Size(); j++) {
  1516. err[i][j] = U_offsurf[i][j] - (U_sl[i*PotentialBasis::Size()+j] + U_dl[i*PotentialBasis::Size()+j]);
  1517. max_err = std::max<Real>(max_err, fabs(err[i][j]));
  1518. }
  1519. }
  1520. { // Print error
  1521. Real glb_err;
  1522. comm.Allreduce(Ptr2ConstItr<Real>(&max_err,1), Ptr2Itr<Real>(&glb_err,1), 1, Comm::CommOp::MAX);
  1523. if (!comm.Rank()) std::cout<<"Error = "<<glb_err<<'\n';
  1524. }
  1525. { // Write VTK output
  1526. VTUData vtu;
  1527. vtu.AddElems(elements_trg, err, ORDER);
  1528. vtu.WriteVTK("err", comm);
  1529. }
  1530. { // Write VTK output
  1531. VTUData vtu;
  1532. vtu.AddElems(elements_trg, U_offsurf, ORDER);
  1533. vtu.WriteVTK("U", comm);
  1534. }
  1535. }
  1536. Profile::print(&comm);
  1537. }
  1538. private:
  1539. static void scan(Vector<Long>& dsp, const Vector<Long>& cnt) {
  1540. dsp.ReInit(cnt.Dim());
  1541. if (cnt.Dim()) dsp[0] = 0;
  1542. omp_par::scan(cnt.begin(), dsp.begin(), cnt.Dim());
  1543. }
  1544. template <class Basis> static void eval_basis(Vector<Real>& value, const Vector<Basis> X, Integer dof, Integer Nnds, const typename Basis::EvalOpType& EvalOp) {
  1545. Long Nelem = X.Dim() / dof;
  1546. SCTL_ASSERT(X.Dim() == Nelem * dof);
  1547. value.ReInit(Nelem*Nnds*dof);
  1548. Matrix<Real> X_(Nelem*dof, Nnds, value.begin(),false);
  1549. Basis::Eval(X_, X, EvalOp);
  1550. for (Long j = 0; j < Nelem; j++) { // Rearrange data
  1551. Matrix<Real> X(Nnds, dof, X_[j*dof], false);
  1552. X = Matrix<Real>(dof, Nnds, X_[j*dof], false).Transpose();
  1553. }
  1554. }
  1555. template <int CoordDim, int FnDim, class FnBasis, class ElemList> static void DiscretizeSurfaceFn(Vector<FnBasis>& U, const ElemList& elements, std::function<void(Real*,Real*,Real*)> fn) {
  1556. using CoordBasis = typename ElemList::CoordBasis;
  1557. const long Nelem = elements.NElem();
  1558. U.ReInit(Nelem * FnDim);
  1559. Matrix<Real> X, X_grad;
  1560. { // Set X, X_grad
  1561. Vector<CoordBasis> coord = elements.ElemVector();
  1562. Vector<CoordBasis> coord_grad;
  1563. CoordBasis::Grad(coord_grad, coord);
  1564. const auto Meval = CoordBasis::SetupEval(FnBasis::Nodes());
  1565. CoordBasis::Eval(X, coord, Meval);
  1566. CoordBasis::Eval(X_grad, coord_grad, Meval);
  1567. }
  1568. for (long i = 0; i < Nelem; i++) {
  1569. for (long j = 0; j < FnBasis::Size(); j++) {
  1570. Real X_[CoordDim], Xn[CoordDim], U_[FnDim];
  1571. for (long k = 0; k < CoordDim; k++) {
  1572. X_[k] = X[i*CoordDim+k][j];
  1573. }
  1574. { // Set Xn
  1575. Real Xu[CoordDim], Xv[CoordDim];
  1576. for (long k = 0; k < CoordDim; k++) {
  1577. Xu[k] = X_grad[(i*CoordDim+k)*2+0][j];
  1578. Xv[k] = X_grad[(i*CoordDim+k)*2+1][j];
  1579. }
  1580. Real dA = 0;
  1581. for (long k = 0; k < CoordDim; k++) {
  1582. Xn[k] = Xu[(k+1)%CoordDim] * Xv[(k+2)%CoordDim];
  1583. Xn[k] -= Xv[(k+1)%CoordDim] * Xu[(k+2)%CoordDim];
  1584. dA += Xn[k] * Xn[k];
  1585. }
  1586. dA = sqrt(dA);
  1587. for (long k = 0; k < CoordDim; k++) {
  1588. Xn[k] /= dA;
  1589. }
  1590. }
  1591. fn(U_, X_, Xn);
  1592. for (long k = 0; k < FnDim; k++) {
  1593. U[i*FnDim+k][j] = U_[k];
  1594. }
  1595. }
  1596. }
  1597. }
  1598. Vector<Real> Xt_;
  1599. Matrix<Real> M_singular;
  1600. Matrix<Real> M_near_singular;
  1601. Vector<Pair<Long,Long>> pair_lst;
  1602. Integer order_direct_;
  1603. Real period_length_;
  1604. Comm comm_;
  1605. };
  1606. } // end namespace
  1607. #endif //_SCTL_BOUNDARY_QUADRATURE_HPP_