boundary_quadrature.hpp 69 KB

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