boundary_quadrature.hpp 70 KB

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