surface_op.txx 34 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842
  1. #include <sctl.hpp>
  2. namespace biest {
  3. template <class Real> SurfaceOp<Real>::SurfaceOp(const sctl::Comm& comm, sctl::Long Nt, sctl::Long Np) {
  4. Init(comm, Nt, Np);
  5. }
  6. template <class Real> SurfaceOp<Real>::SurfaceOp(const SurfaceOp& Op) {
  7. Init(Op.comm_, Op.Nt_, Op.Np_);
  8. }
  9. template <class Real> SurfaceOp<Real>& SurfaceOp<Real>::operator=(const SurfaceOp& Op) {
  10. Init(Op.comm_, Op.Nt_, Op.Np_);
  11. return *this;
  12. }
  13. template <class Real> void SurfaceOp<Real>::Upsample(const sctl::Vector<Real>& X0_, sctl::Long Nt0, sctl::Long Np0, sctl::Vector<Real>& X1_, sctl::Long Nt1, sctl::Long Np1) {
  14. auto FFT_Helper = [](const sctl::FFT<Real>& fft, const sctl::Vector<Real>& in, sctl::Vector<Real>& out) {
  15. sctl::Long dof = in.Dim() / fft.Dim(0);
  16. if (out.Dim() != dof * fft.Dim(1)) {
  17. out.ReInit(dof * fft.Dim(1));
  18. }
  19. for (sctl::Long i = 0; i < dof; i++) {
  20. const sctl::Vector<Real> in_(fft.Dim(0), (sctl::Iterator<Real>)in.begin() + i * fft.Dim(0), false);
  21. sctl::Vector<Real> out_(fft.Dim(1), out.begin() + i * fft.Dim(1), false);
  22. fft.Execute(in_, out_);
  23. }
  24. };
  25. assert(X0_.Dim() % (Nt0 * Np0) == 0);
  26. sctl::Long Nt0_ = Nt0;
  27. sctl::Long Np0_ = Np0 / 2 + 1;
  28. sctl::Long Nt1_ = Nt1;
  29. sctl::Long Np1_ = Np1 / 2 + 1;
  30. sctl::FFT<Real> fft_r2c0, fft_c2r_;
  31. { // Initialize fft_r2c0, fft_c2r_
  32. sctl::StaticArray<sctl::Long, 2> fft_dim0 = {Nt0, Np0};
  33. sctl::StaticArray<sctl::Long, 2> fft_dim_ = {Nt1, Np1};
  34. fft_r2c0.Setup(sctl::FFT_Type::R2C, 1, sctl::Vector<sctl::Long>(2, fft_dim0, false), omp_get_max_threads());
  35. fft_c2r_.Setup(sctl::FFT_Type::C2R, 1, sctl::Vector<sctl::Long>(2, fft_dim_, false), omp_get_max_threads());
  36. }
  37. sctl::Vector<Real> tmp0, tmp_;
  38. auto X0__ = X0_; // TODO: make unaligned plans and remove this workaround
  39. FFT_Helper(fft_r2c0, X0__, tmp0);
  40. sctl::Long dof = tmp0.Dim() / (Nt0_*Np0_*2);
  41. SCTL_ASSERT(tmp0.Dim() == dof * Nt0_ * Np0_ * 2);
  42. if (tmp_.Dim() != dof * Nt1_ * Np1_ * 2) tmp_.ReInit(dof * Nt1_ * Np1_ * 2);
  43. tmp_.SetZero();
  44. Real scale = sctl::sqrt<Real>(Nt1 * Np1) / sctl::sqrt<Real>(Nt0 * Np0);
  45. sctl::Long Ntt = std::min(Nt0_, Nt1_);
  46. sctl::Long Npp = std::min(Np0_, Np1_);
  47. for (sctl::Long k = 0; k < dof; k++) {
  48. for (sctl::Long t = 0; t < (Ntt + 1) / 2; t++) {
  49. for (sctl::Long p = 0; p < Npp; p++) {
  50. tmp_[((k * Nt1_ + t) * Np1_ + p) * 2 + 0] = scale * tmp0[((k * Nt0_ + t) * Np0_ + p) * 2 + 0];
  51. tmp_[((k * Nt1_ + t) * Np1_ + p) * 2 + 1] = scale * tmp0[((k * Nt0_ + t) * Np0_ + p) * 2 + 1];
  52. }
  53. }
  54. for (sctl::Long t = 0; t < Ntt / 2; t++) {
  55. for (sctl::Long p = 0; p < Npp; p++) {
  56. tmp_[((k * Nt1_ + (Nt1_ - t - 1)) * Np1_ + p) * 2 + 0] = scale * tmp0[((k * Nt0_ + (Nt0_ - t - 1)) * Np0_ + p) * 2 + 0];
  57. tmp_[((k * Nt1_ + (Nt1_ - t - 1)) * Np1_ + p) * 2 + 1] = scale * tmp0[((k * Nt0_ + (Nt0_ - t - 1)) * Np0_ + p) * 2 + 1];
  58. }
  59. }
  60. }
  61. if (X1_.Dim() != dof * Nt1 * Np1) X1_.ReInit(dof * Nt1 * Np1);
  62. FFT_Helper(fft_c2r_, tmp_, X1_);
  63. { // Floating-point correction
  64. sctl::Long Ut = Nt1 / Nt0;
  65. sctl::Long Up = Np1 / Np0;
  66. if (Nt1 == Nt0 * Ut && Np1 == Np0 * Up) {
  67. for (sctl::Long k = 0; k < dof; k++) {
  68. for (sctl::Long t = 0; t < Nt0; t++) {
  69. for (sctl::Long p = 0; p < Np0; p++) {
  70. X1_[(k * Nt1 + t * Ut) * Np1 + p * Up] = X0_[(k * Nt0 + t) * Np0 + p];
  71. }
  72. }
  73. }
  74. }
  75. }
  76. }
  77. template <class Real> void SurfaceOp<Real>::Grad2D(sctl::Vector<Real>& dX, const sctl::Vector<Real>& X) const {
  78. sctl::Long dof = X.Dim() / (Nt_ * Np_);
  79. assert(X.Dim() == dof * Nt_ * Np_);
  80. if (dX.Dim() != dof * 2 * Nt_ * Np_) {
  81. dX.ReInit(dof * 2 * Nt_ * Np_);
  82. }
  83. sctl::Long Nt = Nt_;
  84. sctl::Long Np = fft_r2c.Dim(1) / (Nt * 2);
  85. SCTL_ASSERT(fft_r2c.Dim(1) == Nt * Np * 2);
  86. sctl::Vector<Real> coeff(fft_r2c.Dim(1));
  87. sctl::Vector<Real> grad_coeff(fft_r2c.Dim(1));
  88. for (sctl::Long k = 0; k < dof; k++) {
  89. fft_r2c.Execute(sctl::Vector<Real>(Nt_*Np_, (sctl::Iterator<Real>)X.begin() + k*Nt_*Np_, false), coeff);
  90. Real scal = (2 * sctl::const_pi<Real>());
  91. #pragma omp parallel for schedule(static)
  92. for (sctl::Long t = 0; t < Nt; t++) { // grad_coeff(t,p) <-- imag * t * coeff(t,p)
  93. for (sctl::Long p = 0; p < Np; p++) {
  94. Real real = coeff[(t * Np + p) * 2 + 0] * scal;
  95. Real imag = coeff[(t * Np + p) * 2 + 1] * scal;
  96. grad_coeff[(t * Np + p) * 2 + 0] = imag * (t - (t > Nt / 2 ? Nt : 0));
  97. grad_coeff[(t * Np + p) * 2 + 1] = -real * (t - (t > Nt / 2 ? Nt : 0));
  98. }
  99. }
  100. { // dX <-- IFFT(grad_coeff)
  101. sctl::Vector<Real> fft_out(Nt_*Np_, dX.begin() + (k*2+0)*Nt_*Np_, false);
  102. fft_c2r.Execute(grad_coeff, fft_out);
  103. }
  104. scal = (2 * sctl::const_pi<Real>());
  105. #pragma omp parallel for schedule(static)
  106. for (sctl::Long t = 0; t < Nt; t++) { // grad_coeff(t,p) <-- imag * p * coeff(t,p)
  107. for (sctl::Long p = 0; p < Np; p++) {
  108. Real real = coeff[(t * Np + p) * 2 + 0] * scal;
  109. Real imag = coeff[(t * Np + p) * 2 + 1] * scal;
  110. grad_coeff[(t * Np + p) * 2 + 0] = imag * p;
  111. grad_coeff[(t * Np + p) * 2 + 1] = -real * p;
  112. }
  113. }
  114. { // dX <-- IFFT(grad_coeff)
  115. sctl::Vector<Real> fft_out(Nt_*Np_, dX.begin() + (k*2+1)*Nt_*Np_, false);
  116. fft_c2r.Execute(grad_coeff, fft_out);
  117. }
  118. }
  119. }
  120. template <class Real> Real SurfaceOp<Real>::SurfNormalAreaElem(sctl::Vector<Real>* normal, sctl::Vector<Real>* area_elem, const sctl::Vector<Real>& dX, const sctl::Vector<Real>* X) const {
  121. if (normal != nullptr && normal->Dim() != COORD_DIM * Nt_ * Np_) {
  122. normal->ReInit(COORD_DIM * Nt_ * Np_);
  123. }
  124. if (area_elem != nullptr && area_elem->Dim() != Nt_ * Np_) {
  125. area_elem->ReInit(Nt_ * Np_);
  126. }
  127. sctl::Long N = Nt_ * Np_;
  128. Real scal = 1 / (Real)N;
  129. if (normal != nullptr && area_elem != nullptr) {
  130. #pragma omp parallel for schedule(static)
  131. for (sctl::Long i = 0; i < N; i++) {
  132. sctl::StaticArray<Real, COORD_DIM> xt = {dX[0*N+i], dX[2*N+i], dX[4*N+i]};
  133. sctl::StaticArray<Real, COORD_DIM> xp = {dX[1*N+i], dX[3*N+i], dX[5*N+i]};
  134. sctl::StaticArray<Real, COORD_DIM> xn;
  135. (*area_elem)[i] = compute_area_elem(xn, xt, xp) * scal;
  136. (*normal)[0*N+i] = xn[0];
  137. (*normal)[1*N+i] = xn[1];
  138. (*normal)[2*N+i] = xn[2];
  139. }
  140. } else if (normal != nullptr) {
  141. #pragma omp parallel for schedule(static)
  142. for (sctl::Long i = 0; i < N; i++) {
  143. sctl::StaticArray<Real, COORD_DIM> xt = {dX[0*N+i], dX[2*N+i], dX[4*N+i]};
  144. sctl::StaticArray<Real, COORD_DIM> xp = {dX[1*N+i], dX[3*N+i], dX[5*N+i]};
  145. sctl::StaticArray<Real, COORD_DIM> xn;
  146. compute_area_elem(xn, xt, xp);
  147. (*normal)[0*N+i] = xn[0];
  148. (*normal)[1*N+i] = xn[1];
  149. (*normal)[2*N+i] = xn[2];
  150. }
  151. } else if (area_elem != nullptr) {
  152. #pragma omp parallel for schedule(static)
  153. for (sctl::Long i = 0; i < N; i++) {
  154. sctl::StaticArray<Real, COORD_DIM> xt = {dX[0*N+i], dX[2*N+i], dX[4*N+i]};
  155. sctl::StaticArray<Real, COORD_DIM> xp = {dX[1*N+i], dX[3*N+i], dX[5*N+i]};
  156. sctl::StaticArray<Real, COORD_DIM> xn;
  157. (*area_elem)[i] = compute_area_elem(xn, xt, xp) * scal;
  158. }
  159. }
  160. Real normal_scal = 0;
  161. if (X != nullptr && normal != nullptr) { // Set Orientation
  162. sctl::Long idx = 0;
  163. for (sctl::Long j = 0; j < Nt_*Np_; j++) {
  164. if ((*X)[idx] < (*X)[j]) idx = j;
  165. }
  166. if ((*normal)[idx] < 0) {
  167. normal_scal = -1;
  168. (*normal) = (*normal) * normal_scal;
  169. } else {
  170. normal_scal = 1;
  171. }
  172. }
  173. return normal_scal;
  174. }
  175. template <class Real> void SurfaceOp<Real>::SurfCurl(sctl::Vector<Real>& CurlF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& normal, const sctl::Vector<Real>& F) const {
  176. sctl::Vector<Real> GradF;
  177. SurfGrad(GradF, dX, F);
  178. { // Compute CurlF <--curl(GradF)
  179. sctl::Long N = Nt_ * Np_;
  180. sctl::Long dof = F.Dim() / (COORD_DIM * Nt_ * Np_);
  181. assert(F.Dim() == dof * COORD_DIM * Nt_ * Np_);
  182. assert(normal.Dim() == COORD_DIM * Nt_ * Np_);
  183. if (CurlF.Dim() != dof * Nt_ * Np_) {
  184. CurlF.ReInit(dof * Nt_ * Np_);
  185. }
  186. #pragma omp parallel for schedule(static)
  187. for (sctl::Long i = 0; i < Nt_ * Np_; i++) {
  188. for (sctl::Long k = 0; k < dof; k++) {
  189. Real curl = 0;
  190. sctl::Long idx = k*COORD_DIM*COORD_DIM + i;
  191. curl += normal[0*N+i] * (GradF[(1*COORD_DIM+2) * N + idx] - GradF[(2*COORD_DIM+1) * N + idx]);
  192. curl += normal[1*N+i] * (GradF[(2*COORD_DIM+0) * N + idx] - GradF[(0*COORD_DIM+2) * N + idx]);
  193. curl += normal[2*N+i] * (GradF[(0*COORD_DIM+1) * N + idx] - GradF[(1*COORD_DIM+0) * N + idx]);
  194. CurlF[k * N + i] = curl;
  195. }
  196. }
  197. }
  198. }
  199. template <class Real> void SurfaceOp<Real>::SurfGrad(sctl::Vector<Real>& GradFvec, const sctl::Vector<Real>& dXvec_, const sctl::Vector<Real>& Fvec) const {
  200. auto transpose0 = [this](sctl::Vector<Real>& V) {
  201. sctl::Long dof = V.Dim() / (Nt_ * Np_);
  202. assert(V.Dim() == dof * Nt_ * Np_);
  203. sctl::Matrix<Real>(Nt_ * Np_, dof, V.begin(), false) = sctl::Matrix<Real>(dof, Nt_ * Np_, V.begin(), false).Transpose();
  204. };
  205. auto transpose1 = [this](sctl::Vector<Real>& V) {
  206. sctl::Long dof = V.Dim() / (Nt_ * Np_);
  207. assert(V.Dim() == Nt_ * Np_ * dof);
  208. sctl::Matrix<Real>(dof, Nt_ * Np_, V.begin(), false) = sctl::Matrix<Real>(Nt_ * Np_, dof, V.begin(), false).Transpose();
  209. };
  210. auto inv2x2 = [](Mat<Real, 2, 2> M) {
  211. Mat<Real, 2, 2> Mout;
  212. Real oodet = 1 / (M[0][0] * M[1][1] - M[0][1] * M[1][0]);
  213. Mout[0][0] = M[1][1] * oodet;
  214. Mout[0][1] = -M[0][1] * oodet;
  215. Mout[1][0] = -M[1][0] * oodet;
  216. Mout[1][1] = M[0][0] * oodet;
  217. return Mout;
  218. };
  219. SCTL_ASSERT(dXvec_.Dim() == 2 * COORD_DIM * Nt_ * Np_);
  220. auto dXvec = dXvec_;
  221. sctl::Vector<Real> dFvec;
  222. Grad2D(dFvec, Fvec);
  223. transpose0(dXvec);
  224. transpose0(dFvec);
  225. { // Set GradFvec
  226. sctl::Long dof = Fvec.Dim() / (Nt_ * Np_);
  227. assert(Fvec.Dim() == Nt_ * Np_ * dof);
  228. if (GradFvec.Dim() != dof * COORD_DIM * Nt_ * Np_) {
  229. GradFvec.ReInit(dof * COORD_DIM * Nt_ * Np_);
  230. }
  231. #pragma omp parallel for schedule(static)
  232. for (sctl::Long i = 0; i < Nt_ * Np_; i++) {
  233. const Mat<Real, 3, 2, false> dX_du(dXvec.begin() + i*3*2);
  234. const auto du_dX = inv2x2(dX_du.Transpose() * dX_du) * dX_du.Transpose();
  235. for (sctl::Long k = 0; k < dof; k++) {
  236. const Mat<Real, 1, 2, false> dF_du(dFvec.begin() + (i*dof+k)*1*2);
  237. Mat<Real, 1, 3, false> GradF(GradFvec.begin() + (i*dof+k)*1*3);
  238. GradF = dF_du * du_dX;
  239. }
  240. }
  241. transpose1(GradFvec);
  242. }
  243. }
  244. template <class Real> void SurfaceOp<Real>::SurfDiv(sctl::Vector<Real>& DivF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& F) const {
  245. sctl::Vector<Real> GradF;
  246. SurfGrad(GradF, dX, F);
  247. { // Compute DivF <-- trace(GradF)
  248. sctl::Long N = Nt_ * Np_;
  249. sctl::Long dof = F.Dim() / (COORD_DIM * Nt_ * Np_);
  250. assert(F.Dim() == dof * COORD_DIM * Nt_ * Np_);
  251. if (DivF.Dim() != dof * Nt_ * Np_) {
  252. DivF.ReInit(dof * Nt_ * Np_);
  253. }
  254. #pragma omp parallel for schedule(static)
  255. for (sctl::Long i = 0; i < Nt_ * Np_; i++) {
  256. for (sctl::Long k = 0; k < dof; k++) {
  257. Real trace = 0;
  258. for (sctl::Long l = 0; l < COORD_DIM; l++) {
  259. trace += GradF[(k * COORD_DIM * COORD_DIM + (COORD_DIM + 1) * l) * N + i];
  260. }
  261. DivF[k * N + i] = trace;
  262. }
  263. }
  264. }
  265. }
  266. template <class Real> void SurfaceOp<Real>::SurfLap(sctl::Vector<Real>& LapF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& F) const {
  267. sctl::Vector<Real> GradF;
  268. SurfGrad(GradF, dX, F);
  269. SurfDiv(LapF, dX, GradF);
  270. }
  271. template <class Real> void SurfaceOp<Real>::SurfInteg(sctl::Vector<Real>& I, const sctl::Vector<Real>& area_elem, const sctl::Vector<Real>& F) const {
  272. sctl::Long dof = F.Dim() / (Nt_ * Np_);
  273. SCTL_ASSERT(F.Dim() == dof * Nt_ * Np_);
  274. SCTL_ASSERT(area_elem.Dim() == Nt_ * Np_);
  275. if (I.Dim() != dof) I.ReInit(dof);
  276. for (sctl::Long k = 0; k < dof; k++) {
  277. Real sum = 0;
  278. #pragma omp parallel for schedule(static) reduction(+:sum)
  279. for (sctl::Long i = 0; i < Nt_ * Np_; i++) {
  280. sum += area_elem[i] * F[k * Nt_ * Np_ + i];
  281. }
  282. I[k] = sum;
  283. }
  284. }
  285. template <class Real> void SurfaceOp<Real>::ProjZeroMean(sctl::Vector<Real>& Fproj, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& F) const {
  286. sctl::Long N = this->Nt_ * this->Np_;
  287. sctl::Vector<Real> area_elem;
  288. SurfNormalAreaElem(nullptr, &area_elem, dX, nullptr);
  289. sctl::Vector<Real> Favg, SurfArea, one(N);
  290. one = 1;
  291. SurfInteg(SurfArea, area_elem, one);
  292. SurfInteg(Favg, area_elem, F);
  293. Favg /= SurfArea[0];
  294. sctl::Long dof = Favg.Dim();
  295. Fproj.ReInit(dof * N);
  296. for (sctl::Long k = 0; k < dof; k++) {
  297. #pragma omp parallel for schedule(static)
  298. for (sctl::Long i = 0; i < N; i++) {
  299. Fproj[k * N + i] = F[k * N + i] - Favg[k];
  300. }
  301. }
  302. }
  303. template <class Real> void SurfaceOp<Real>::InvSurfLap(sctl::Vector<Real>& InvLapF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& d2X, const sctl::Vector<Real>& F, Real tol, sctl::Integer max_iter, Real upsample) const {
  304. auto spectral_InvSurfLap = [this](sctl::Vector<Real>& InvLapF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& F, sctl::Long Nt_, sctl::Long Np_, Real tol, sctl::Integer max_iter) {
  305. sctl::FFT<Real> fft_r2c, fft_c2r;
  306. sctl::StaticArray<sctl::Long, 2> fft_dim = {Nt_, Np_};
  307. fft_r2c.Setup(sctl::FFT_Type::R2C, 1, sctl::Vector<sctl::Long>(2, fft_dim, false), omp_get_max_threads());
  308. fft_c2r.Setup(sctl::FFT_Type::C2R, 1, sctl::Vector<sctl::Long>(2, fft_dim, false), omp_get_max_threads());
  309. std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> IdentityOp =[](sctl::Vector<Real>& Sx, const sctl::Vector<Real>& x) { Sx = x; };
  310. std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> FPrecond = [Nt_, Np_, &fft_r2c, &fft_c2r, &dX](sctl::Vector<Real>& Sx, const sctl::Vector<Real>& x) {
  311. sctl::Vector<Real> coeff;
  312. fft_r2c.Execute(x, coeff);
  313. Real Lt = 0, Lp = 0;
  314. sctl::Long N = Nt_ * Np_;
  315. for (sctl::Long i = 0; i < N; i++) {
  316. Real Xt = dX[(2*0+0) * N + i];
  317. Real Yt = dX[(2*1+0) * N + i];
  318. Real Zt = dX[(2*2+0) * N + i];
  319. Real Xp = dX[(2*0+1) * N + i];
  320. Real Yp = dX[(2*1+1) * N + i];
  321. Real Zp = dX[(2*2+1) * N + i];
  322. Lt += sqrt(Xt * Xt + Yt * Yt + Zt * Zt);
  323. Lp += sqrt(Xp * Xp + Yp * Yp + Zp * Zp);
  324. }
  325. Real invLt2 = N * N / (Lt * Lt);
  326. Real invLp2 = N * N / (Lp * Lp);
  327. sctl::Long Nt0_ = Nt_;
  328. sctl::Long Np0_ = fft_r2c.Dim(1) / (Nt0_ * 2);
  329. for (sctl::Long t = 0; t < Nt0_; t++) { // coeff(t,p) <-- coeff(t,p) / (t*t + p*p)
  330. for (sctl::Long p = 0; p < Np0_; p++) {
  331. sctl::Long tt = (t - (t > Nt0_ / 2 ? Nt0_ : 0));
  332. sctl::Long pp = p;
  333. Real scal = (tt || pp ? 1 / (Real)(tt*tt*invLt2 + pp*pp*invLp2) : 0);
  334. coeff[(t * Np0_ + p) * 2 + 0] *= scal;
  335. coeff[(t * Np0_ + p) * 2 + 1] *= scal;
  336. }
  337. }
  338. fft_c2r.Execute(coeff, Sx);
  339. };
  340. SurfaceOp<Real> Sop(comm_, Nt_, Np_);
  341. Sop.InvSurfLapPrecond(InvLapF, FPrecond, IdentityOp, dX, F, tol, max_iter);
  342. };
  343. sctl::Long Nt_up = (sctl::Long)(Nt_ * upsample);
  344. sctl::Long Np_up = (sctl::Long)(Np_ * upsample);
  345. sctl::Vector<Real> InvLapF_up, dX_up, F_up;
  346. Upsample(F, Nt_, Np_, F_up, Nt_up, Np_up);
  347. Upsample(dX, Nt_, Np_, dX_up, Nt_up, Np_up);
  348. spectral_InvSurfLap(InvLapF_up, dX_up, F_up, Nt_up, Np_up, tol, max_iter);
  349. Upsample(InvLapF_up, Nt_up, Np_up, InvLapF, Nt_, Np_);
  350. }
  351. template <class Real> void SurfaceOp<Real>::GradInvSurfLap(sctl::Vector<Real>& GradInvLapF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& d2X, const sctl::Vector<Real>& F, Real tol, sctl::Integer max_iter, Real upsample) const {
  352. auto spectral_GradInvSurfLap = [this](sctl::Vector<Real>& GradInvLapF, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& F, sctl::Long Nt_, sctl::Long Np_, Real tol, sctl::Integer max_iter) {
  353. sctl::FFT<Real> fft_r2c, fft_c2r;
  354. sctl::StaticArray<sctl::Long, 2> fft_dim = {Nt_, Np_};
  355. fft_r2c.Setup(sctl::FFT_Type::R2C, 1, sctl::Vector<sctl::Long>(2, fft_dim, false), omp_get_max_threads());
  356. fft_c2r.Setup(sctl::FFT_Type::C2R, 1, sctl::Vector<sctl::Long>(2, fft_dim, false), omp_get_max_threads());
  357. std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> IdentityOp = [](sctl::Vector<Real>& Sx, const sctl::Vector<Real>& x) { Sx = x; };
  358. std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> FPrecond = [Nt_, Np_, &fft_r2c, &fft_c2r, &dX](sctl::Vector<Real>& Sx, const sctl::Vector<Real>& x) {
  359. sctl::Vector<Real> coeff;
  360. fft_r2c.Execute(x, coeff);
  361. Real Lt = 0, Lp = 0;
  362. sctl::Long N = Nt_ * Np_;
  363. #pragma omp parallel for schedule(static) reduction(+:Lt,Lp)
  364. for (sctl::Long i = 0; i < N; i++) {
  365. Real Xt = dX[(2*0+0) * N + i];
  366. Real Yt = dX[(2*1+0) * N + i];
  367. Real Zt = dX[(2*2+0) * N + i];
  368. Real Xp = dX[(2*0+1) * N + i];
  369. Real Yp = dX[(2*1+1) * N + i];
  370. Real Zp = dX[(2*2+1) * N + i];
  371. Lt += sqrt(Xt * Xt + Yt * Yt + Zt * Zt);
  372. Lp += sqrt(Xp * Xp + Yp * Yp + Zp * Zp);
  373. }
  374. Real invLt2 = N * N / (Lt * Lt);
  375. Real invLp2 = N * N / (Lp * Lp);
  376. sctl::Long Nt0_ = Nt_;
  377. sctl::Long Np0_ = fft_r2c.Dim(1) / (Nt0_ * 2);
  378. #pragma omp parallel for schedule(static)
  379. for (sctl::Long t = 0; t < Nt0_; t++) { // coeff(t,p) <-- coeff(t,p) / (t*t + p*p)
  380. for (sctl::Long p = 0; p < Np0_; p++) {
  381. sctl::Long tt = (t - (t > Nt0_ / 2 ? Nt0_ : 0));
  382. sctl::Long pp = p;
  383. Real scal = (tt || pp ? 1 / (Real)(tt*tt*invLt2 + pp*pp*invLp2) : 0);
  384. coeff[(t * Np0_ + p) * 2 + 0] *= scal;
  385. coeff[(t * Np0_ + p) * 2 + 1] *= scal;
  386. }
  387. }
  388. fft_c2r.Execute(coeff, Sx);
  389. };
  390. sctl::Vector<Real> InvLapF;
  391. SurfaceOp<Real> Sop(comm_, Nt_, Np_);
  392. Sop.InvSurfLapPrecond(InvLapF, FPrecond, IdentityOp, dX, F, tol, max_iter);
  393. Sop.SurfGrad(GradInvLapF, dX, InvLapF);
  394. };
  395. sctl::Long Nt_up = (sctl::Long)(Nt_ * upsample);
  396. sctl::Long Np_up = (sctl::Long)(Np_ * upsample);
  397. sctl::Vector<Real> GradInvLapF_up, dX_up, F_up;
  398. Upsample(F, Nt_, Np_, F_up, Nt_up, Np_up);
  399. Upsample(dX, Nt_, Np_, dX_up, Nt_up, Np_up);
  400. spectral_GradInvSurfLap(GradInvLapF_up, dX_up, F_up, Nt_up, Np_up, tol, max_iter);
  401. Upsample(GradInvLapF_up, Nt_up, Np_up, GradInvLapF, Nt_, Np_);
  402. }
  403. template <class Real> void SurfaceOp<Real>::InvSurfLapPrecond(sctl::Vector<Real>& InvLapF, std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> LeftPrecond, std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> RightPrecond, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& F, Real tol, sctl::Integer max_iter) const {
  404. typename sctl::ParallelSolver<Real>::ParallelOp fn = [this,&LeftPrecond,&RightPrecond,&dX](sctl::Vector<Real>* fx, const sctl::Vector<Real>& x) {
  405. (*fx) = 0;
  406. SCTL_ASSERT(fx);
  407. sctl::Vector<Real> Sx, LSx;
  408. RightPrecond(Sx, x);
  409. SurfLap(LSx, dX, Sx);
  410. sctl::Vector<Real> Proj_Sx;
  411. ProjZeroMean(Proj_Sx, dX, Sx);
  412. LSx += (Sx - Proj_Sx)*1;
  413. LeftPrecond(*fx, LSx);
  414. };
  415. sctl::Vector<Real> Fproj, SFproj, InvSInvLapF, InvLapF_;
  416. ProjZeroMean(Fproj, dX, F);
  417. LeftPrecond(SFproj, Fproj);
  418. solver(&InvSInvLapF, fn, SFproj, tol, max_iter);
  419. RightPrecond(InvLapF_, InvSInvLapF);
  420. ProjZeroMean(InvLapF, dX, InvLapF_);
  421. }
  422. template <class Real> template <sctl::Integer KDIM0, sctl::Integer KDIM1> void SurfaceOp<Real>::EvalSurfInteg(sctl::Vector<Real>& Utrg, const sctl::Vector<Real>& Xtrg, const sctl::Vector<Real>& Xsrc, const sctl::Vector<Real>& Xn_src, const sctl::Vector<Real>& Xa_src, const sctl::Vector<Real>& Fsrc, const KernelFunction<Real,COORD_DIM,KDIM0,KDIM1>& ker) const {
  423. sctl::Long Nsrc = Nt_* Np_;
  424. sctl::Long dof = Fsrc.Dim() / (KDIM0 * Nsrc);
  425. SCTL_ASSERT( Xsrc.Dim() == COORD_DIM * Nsrc);
  426. SCTL_ASSERT(Xn_src.Dim() == COORD_DIM * Nsrc);
  427. SCTL_ASSERT(Xa_src.Dim() == Nsrc);
  428. SCTL_ASSERT( Fsrc.Dim() == dof*KDIM0 * Nsrc);
  429. sctl::Vector<Real> Fa(dof * KDIM0 * Nsrc);
  430. for (sctl::Long i = 0; i < Nsrc; i++) {
  431. for (sctl::Integer j = 0; j < dof * KDIM0; j++) {
  432. Fa[j * Nsrc + i] = Fsrc[j * Nsrc + i] * Xa_src[i];
  433. }
  434. }
  435. sctl::Long np = comm_.Size();
  436. sctl::Long rank = comm_.Rank();
  437. sctl::Long Ntrg = Xtrg.Dim() / COORD_DIM;
  438. SCTL_ASSERT(Xtrg.Dim() == COORD_DIM * Ntrg);
  439. sctl::Long a = (rank + 0) * Ntrg / np;
  440. sctl::Long b = (rank + 1) * Ntrg / np;
  441. sctl::Vector<Real> TrgX, TrgU;
  442. { // Init TrgX, TrgU
  443. TrgX.ReInit(COORD_DIM * (b - a));
  444. TrgU.ReInit(dof * KDIM1 * (b - a));
  445. for (sctl::Integer k = 0; k < COORD_DIM; k++) {
  446. for (sctl::Long i = a; i < b; i++) {
  447. TrgX[k * (b - a) + i - a] = Xtrg[k * Ntrg + i];
  448. }
  449. }
  450. TrgU = 0;
  451. }
  452. ker(Xsrc, Xn_src, Fa, TrgX, TrgU);
  453. { // Set Utrg <-- Reduce(TrgU)
  454. sctl::Vector<Real> Uloc(dof * KDIM1 * Ntrg), Uglb(dof * KDIM1 * Ntrg); Uloc = 0;
  455. for (sctl::Integer k = 0; k < dof * ker.Dim(1); k++) {
  456. for (sctl::Long i = a; i < b; i++) {
  457. Uloc[k * Ntrg + i] = TrgU[k * (b - a) + i - a];
  458. }
  459. }
  460. comm_.template Allreduce<Real>(Uloc.begin(), Uglb.begin(), Uglb.Dim(), sctl::Comm::CommOp::SUM);
  461. if (Utrg.Dim() != dof * KDIM1 * Ntrg) {
  462. Utrg.ReInit(dof * KDIM1 * Ntrg);
  463. Utrg = 0;
  464. }
  465. Utrg += Uglb;
  466. }
  467. }
  468. template <class Real> template <class SingularCorrection, class Kernel> void SurfaceOp<Real>::SetupSingularCorrection(sctl::Vector<SingularCorrection>& singular_correction, sctl::Integer TRG_SKIP, const sctl::Vector<Real>& Xsrc, const sctl::Vector<Real>& dXsrc, const Kernel& ker, Real normal_scal) const {
  469. sctl::Long np = comm_.Size();
  470. sctl::Long rank = comm_.Rank();
  471. sctl::Long Nt0 = Nt_ / TRG_SKIP;
  472. sctl::Long Np0 = Np_ / TRG_SKIP;
  473. sctl::Long a = (rank + 0) * Nt0 * Np0 / np;
  474. sctl::Long b = (rank + 1) * Nt0 * Np0 / np;
  475. singular_correction.ReInit(b - a);
  476. #pragma omp parallel for schedule(static)
  477. for (sctl::Long i = a; i < b; i++) {
  478. sctl::Long t = (i / Np0) * TRG_SKIP;
  479. sctl::Long p = (i % Np0) * TRG_SKIP;
  480. singular_correction[i - a].Setup(TRG_SKIP, Nt_, Np_, Xsrc, dXsrc, t, p, ker, normal_scal);
  481. }
  482. }
  483. template <class Real> template <class SingularCorrection> void SurfaceOp<Real>::EvalSingularCorrection(sctl::Vector<Real>& U, const sctl::Vector<SingularCorrection>& singular_correction, sctl::Integer kdim0, sctl::Integer kdim1, const sctl::Vector<Real>& F) const {
  484. sctl::Long Ntrg;
  485. { // Set N
  486. sctl::Long Nloc = singular_correction.Dim();
  487. comm_.Allreduce(sctl::Ptr2ConstItr<sctl::Long>(&Nloc,1), sctl::Ptr2Itr<sctl::Long>(&Ntrg,1), 1, sctl::Comm::CommOp::SUM);
  488. }
  489. sctl::Long dof = F.Dim() / (kdim0 * Nt_ * Np_);
  490. SCTL_ASSERT(F.Dim() == dof * kdim0 * Nt_ * Np_);
  491. sctl::Vector<Real> Uloc(dof * kdim1 * Ntrg), Uglb(dof * kdim1 * Ntrg);
  492. for (sctl::Long k = 0; k < dof; k++) { // Set Uloc
  493. const sctl::Vector<Real> F_(kdim0 * Nt_ * Np_, (sctl::Iterator<Real>)F.begin() + k * kdim0 * Nt_ * Np_, false);
  494. sctl::Vector<Real> U_(kdim1 * Ntrg, Uloc.begin() + k * kdim1 * Ntrg, false);
  495. U_ = 0;
  496. #pragma omp parallel for schedule(static)
  497. for (sctl::Long i = 0; i < singular_correction.Dim(); i++) {
  498. singular_correction[i](F_, U_);
  499. }
  500. }
  501. comm_.Allreduce(Uloc.begin(), Uglb.begin(), Uglb.Dim(), sctl::Comm::CommOp::SUM);
  502. if (U.Dim() != dof * kdim1 * Ntrg) { // Init U
  503. U.ReInit(dof * kdim1 * Ntrg);
  504. U = 0;
  505. }
  506. U += Uglb;
  507. }
  508. template <class Real> void SurfaceOp<Real>::HodgeDecomp(sctl::Vector<Real>& Vn, sctl::Vector<Real>& Vd, sctl::Vector<Real>& Vc, sctl::Vector<Real>& Vh, const sctl::Vector<Real>& V, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& d2X, const sctl::Vector<Real>& normal, Real tol, sctl::Long max_iter) const {
  509. sctl::Long N = Nt_ * Np_;
  510. SCTL_ASSERT(V.Dim() == COORD_DIM * N);
  511. if (Vn.Dim() != COORD_DIM * N) Vn.ReInit(COORD_DIM * N);
  512. if (Vd.Dim() != COORD_DIM * N) Vd.ReInit(COORD_DIM * N);
  513. if (Vc.Dim() != COORD_DIM * N) Vc.ReInit(COORD_DIM * N);
  514. if (Vh.Dim() != COORD_DIM * N) Vh.ReInit(COORD_DIM * N);
  515. { // Set Vn
  516. Vn.ReInit(V.Dim());
  517. for (sctl::Long i = 0; i < N; i++) {
  518. Real VdotN = 0;
  519. for (sctl::Long k = 0; k < COORD_DIM; k++) VdotN += V[k * N + i] * normal[k * N + i];
  520. for (sctl::Long k = 0; k < COORD_DIM; k++) Vn[k * N + i] = VdotN * normal[k * N + i];
  521. }
  522. }
  523. { // Set Vd = Grad(InvLap(Div(V)))
  524. sctl::Vector<Real> DivV, GradInvLapDivV;
  525. SurfDiv(DivV, dX, V);
  526. GradInvSurfLap(Vd, dX, d2X, DivV, tol * max_norm(V) / max_norm(DivV), max_iter, 1.5);
  527. }
  528. { // Set Vc = n x Grad(InvLap(Div(V x n)))
  529. auto cross_prod = [](sctl::Vector<Real>& axb, const sctl::Vector<Real>& a, const sctl::Vector<Real>& b) {
  530. sctl::Long N = a.Dim() / COORD_DIM;
  531. SCTL_ASSERT(a.Dim() == COORD_DIM * N);
  532. SCTL_ASSERT(b.Dim() == COORD_DIM * N);
  533. if (axb.Dim() != COORD_DIM * N) axb.ReInit(COORD_DIM * N);
  534. for (sctl::Long i = 0; i < N; i++) {
  535. axb[0 * N + i] = a[1 * N + i] * b[2 * N + i] - a[2 * N + i] * b[1 * N + i];
  536. axb[1 * N + i] = a[2 * N + i] * b[0 * N + i] - a[0 * N + i] * b[2 * N + i];
  537. axb[2 * N + i] = a[0 * N + i] * b[1 * N + i] - a[1 * N + i] * b[0 * N + i];
  538. }
  539. };
  540. sctl::Vector<Real> Vxn, DivVxn, GradInvLapDivVxn;
  541. cross_prod(Vxn, V, normal);
  542. SurfDiv(DivVxn, dX, Vxn);
  543. GradInvSurfLap(GradInvLapDivVxn, dX, d2X, DivVxn, tol * max_norm(V) / max_norm(DivVxn), max_iter, 1.5);
  544. cross_prod(Vc, normal, GradInvLapDivVxn);
  545. }
  546. Vh = V - Vn - Vd - Vc;
  547. }
  548. template <class Real> void SurfaceOp<Real>::test_SurfGrad(sctl::Long Nt, sctl::Long Np, SurfType surf_type, const sctl::Comm& comm) {
  549. Surface<Real> S(Nt, Np, surf_type);
  550. const sctl::Vector<Real> f(Nt * Np, (sctl::Iterator<Real>)S.Coord().begin() + 1 * Nt * Np);
  551. sctl::Vector<Real> Lf, dX, Sn;
  552. SurfaceOp Op(comm, Nt, Np);
  553. Op.Grad2D(dX, S.Coord());
  554. Op.SurfGrad(Lf, dX, f);
  555. Op.SurfNormalAreaElem(&Sn, nullptr, dX, &S.Coord());
  556. for (sctl::Long i = 0; i < Nt; i++) { // Lf <-- Lf - (e1 - e1 . n)
  557. for (sctl::Long j = 0; j < Np; j++) {
  558. Real n[3];
  559. n[0] = Sn[(0 * Nt + i) * Np + j];
  560. n[1] = Sn[(1 * Nt + i) * Np + j];
  561. n[2] = Sn[(2 * Nt + i) * Np + j];
  562. Real x[3] = {0,1,0};
  563. Real ndotx = x[0]*n[0] + x[1]*n[1] + x[2]*n[2];
  564. Lf[(0 * Nt + i) * Np + j] -= x[0] - ndotx * n[0];
  565. Lf[(1 * Nt + i) * Np + j] -= x[1] - ndotx * n[1];
  566. Lf[(2 * Nt + i) * Np + j] -= x[2] - ndotx * n[2];
  567. }
  568. }
  569. std::cout<<"|e|_inf = "<<max_norm(Lf)<<'\n';
  570. }
  571. template <class Real> void SurfaceOp<Real>::test_SurfLap(sctl::Long Nt, sctl::Long Np, SurfType surf_type, const sctl::Comm& comm) {
  572. Surface<Real> S(Nt, Np, surf_type);
  573. sctl::Vector<Real> dX, d2X;
  574. SurfaceOp Op(comm, Nt, Np);
  575. Op.Grad2D(dX, S.Coord());
  576. Op.Grad2D(d2X, dX);
  577. sctl::Vector<Real> u0, f0, f1;
  578. Op.LaplaceBeltramiReference(f0, u0, S.Coord(), dX, d2X);
  579. Op.SurfLap(f1, dX, u0);
  580. //WriteVTK("u0", S, u0);
  581. //WriteVTK("f0", S, f0);
  582. //WriteVTK("f1", S, f1);
  583. //WriteVTK("err", S, f1-f0);
  584. { // print error
  585. std::cout<<"|u|_inf = "<<max_norm(u0)<<'\n';
  586. std::cout<<"|f|_inf = "<<max_norm(f0)<<'\n';
  587. std::cout<<"|e|_inf = "<<max_norm(f1-f0)<<'\n';
  588. }
  589. }
  590. template <class Real> void SurfaceOp<Real>::test_InvSurfLap(sctl::Long Nt, sctl::Long Np, SurfType surf_type, Real gmres_tol, sctl::Long gmres_iter, std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)>* LeftPrecond, std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)>* RightPrecond, const sctl::Comm& comm) {
  591. Surface<Real> S(Nt, Np, surf_type);
  592. sctl::Vector<Real> dX, d2X;
  593. SurfaceOp Op(comm, Nt, Np);
  594. Op.Grad2D(dX, S.Coord());
  595. Op.Grad2D(d2X, dX);
  596. sctl::Vector<Real> f0, u0, u1;
  597. Op.LaplaceBeltramiReference(f0, u0, S.Coord(), dX, d2X);
  598. sctl::Profile::Tic("Solve", &comm);
  599. if (LeftPrecond == nullptr && RightPrecond == nullptr) {
  600. Op.InvSurfLap(u1, dX, d2X, f0, gmres_tol, gmres_iter, 1.0);
  601. } else {
  602. std::function<void(sctl::Vector<Real>&, const sctl::Vector<Real>&)> IdentityOp = [](sctl::Vector<Real>& Sx, const sctl::Vector<Real>& x) { Sx = x; };
  603. auto& LPrecond = (LeftPrecond ? *LeftPrecond : IdentityOp);
  604. auto& RPrecond = (RightPrecond ? *RightPrecond : IdentityOp);
  605. Op.InvSurfLapPrecond(u1, LPrecond, RPrecond, dX, f0, gmres_tol, gmres_iter);
  606. }
  607. sctl::Profile::Toc();
  608. //WriteVTK("u0", S, u0);
  609. //WriteVTK("u1", S, u1);
  610. //WriteVTK("err", S, u1-u0);
  611. //WriteVTK("f", S, f0);
  612. { // print error
  613. auto inf_norm = [](const sctl::Vector<Real>& v){
  614. Real max_err = 0;
  615. for (const auto& x: v) max_err = std::max(max_err, fabs(x));
  616. return max_err;
  617. };
  618. std::cout<<"|f|_inf = "<<inf_norm(f0)<<'\n';
  619. std::cout<<"|u|_inf = "<<inf_norm(u0)<<'\n';
  620. std::cout<<"|e|_inf = "<<inf_norm(u1-u0)<<'\n';
  621. }
  622. }
  623. template <class Real> void SurfaceOp<Real>::test_HodgeDecomp(sctl::Long Nt, sctl::Long Np, SurfType surf_type, const sctl::Comm& comm) {
  624. sctl::Long N = Nt * Np;
  625. Surface<Real> S(Nt, Np, surf_type);
  626. sctl::Vector<Real> dX, d2X, normal, area_elem;
  627. SurfaceOp Op(comm, Nt, Np);
  628. Op.Grad2D(dX, S.Coord());
  629. Op.Grad2D(d2X, dX);
  630. Op.SurfNormalAreaElem(&normal, &area_elem, dX, &S.Coord());
  631. Real tol = 1e-9;
  632. sctl::Long max_iter = 100;
  633. sctl::Vector<Real> V(COORD_DIM * N), Vn, Vd, Vc, Vh;
  634. for (sctl::Long i = 0; i < N; i++) { // Set V = dX/du
  635. for (sctl::Long k = 0; k < COORD_DIM; k++) {
  636. V[k * N + i] = dX[(2*k+0) * N + i];
  637. }
  638. }
  639. Op.HodgeDecomp(Vn, Vd, Vc, Vh, V, dX, d2X, normal, tol, max_iter);
  640. auto l2norm = [](const sctl::Vector<Real>& v) {
  641. Real sum = 0;
  642. for (const auto x : v) sum += x*x;
  643. return sqrt(sum);
  644. };
  645. std::cout<<"|Vn|_2: "<<l2norm(Vn)<<'\n';
  646. std::cout<<"|Vd|_2: "<<l2norm(Vd)<<'\n';
  647. std::cout<<"|Vc|_2: "<<l2norm(Vc)<<'\n';
  648. std::cout<<"|Vh|_2: "<<l2norm(Vh)<<'\n';
  649. std::cout<<"|V |_2: "<<l2norm(V )<<'\n';
  650. //WriteVTK("V" , S, V);
  651. //WriteVTK("Vn", S, Vn);
  652. //WriteVTK("Vd", S, Vd);
  653. //WriteVTK("Vc", S, Vc);
  654. //WriteVTK("Vh", S, Vh);
  655. }
  656. template <class Real> void SurfaceOp<Real>::Init(const sctl::Comm& comm, sctl::Long Nt, sctl::Long Np) {
  657. Nt_ = Nt;
  658. Np_ = Np;
  659. comm_ = comm;
  660. solver = sctl::ParallelSolver<Real>(comm_,false);
  661. if (!Nt_ || !Np_) return;
  662. sctl::StaticArray<sctl::Long, 2> fft_dim = {Nt_, Np_};
  663. fft_r2c.Setup(sctl::FFT_Type::R2C, 1, sctl::Vector<sctl::Long>(2, fft_dim, false));
  664. fft_c2r.Setup(sctl::FFT_Type::C2R, 1, sctl::Vector<sctl::Long>(2, fft_dim, false));
  665. }
  666. template <class Real> Real SurfaceOp<Real>::max_norm(const sctl::Vector<Real>& x) {
  667. Real err = 0;
  668. for (const auto& a : x) err = std::max(err, sctl::fabs<Real>(a));
  669. return err;
  670. }
  671. template <class Real> Real SurfaceOp<Real>::compute_area_elem(sctl::StaticArray<Real, COORD_DIM>& xn, const sctl::StaticArray<Real, COORD_DIM>& xt, const sctl::StaticArray<Real, COORD_DIM>& xp) {
  672. xn[0] = xt[1] * xp[2] - xp[1] * xt[2];
  673. xn[1] = xt[2] * xp[0] - xp[2] * xt[0];
  674. xn[2] = xt[0] * xp[1] - xp[0] * xt[1];
  675. Real xa = sqrt(xn[0] * xn[0] + xn[1] * xn[1] + xn[2] * xn[2]);
  676. Real xa_inv = 1 / xa;
  677. xn[0] *= xa_inv;
  678. xn[1] *= xa_inv;
  679. xn[2] *= xa_inv;
  680. return xa;
  681. };
  682. template <class Real> void SurfaceOp<Real>::LaplaceBeltramiReference(sctl::Vector<Real>& f0, sctl::Vector<Real>& u0, const sctl::Vector<Real>& X, const sctl::Vector<Real>& dX, const sctl::Vector<Real>& d2X) const {
  683. sctl::Vector<Real> normal, area_elem;
  684. SurfNormalAreaElem(&normal, &area_elem, dX, &X);
  685. sctl::StaticArray<Real, COORD_DIM> Xs = {2, 2, 2};
  686. { // Set u0
  687. sctl::Vector<Real> SrcCoord, SrcNormal, SrcValue;
  688. SrcCoord.ReInit(COORD_DIM, Xs, false);;
  689. SrcNormal = SrcCoord;
  690. SrcValue.PushBack(1);
  691. const auto& ker = Laplace3D<Real>::FxU();
  692. ker(SrcCoord, SrcNormal, SrcValue, X, u0);
  693. { // Project u0
  694. sctl::Vector<Real> one, surf_area, u0_avg;
  695. one = u0 * 0 + 1;
  696. SurfInteg(surf_area, area_elem, one);
  697. SurfInteg(u0_avg, area_elem, u0);
  698. u0_avg[0] /= surf_area[0];
  699. u0 -= u0_avg[0];
  700. }
  701. }
  702. { // Set f0
  703. sctl::Long N = Nt_ * Np_;
  704. f0.ReInit(N);
  705. for (sctl::Long i = 0; i < N; i++) {
  706. sctl::StaticArray<Real, COORD_DIM> Xt = {X[0*N+i],X[1*N+i],X[2*N+i]};
  707. sctl::StaticArray<Real, COORD_DIM> Xn = {normal[0*N+i],normal[1*N+i],normal[2*N+i]};
  708. sctl::StaticArray<Real, COORD_DIM> dR = {Xt[0]-Xs[0], Xt[1]-Xs[1], Xt[2]-Xs[2]};
  709. Real R2 = dR[0]*dR[0] + dR[1]*dR[1] + dR[2]*dR[2];
  710. Real invR = 1.0/sqrt(R2);
  711. Real invR3 = invR*invR*invR;
  712. Real invR5 = invR3*invR*invR;
  713. static const Real scal = 1 / (4 * sctl::const_pi<Real>());
  714. sctl::StaticArray<Real, COORD_DIM> dG = {
  715. -scal*dR[0]*invR3, -scal*dR[1]*invR3, -scal*dR[2]*invR3
  716. };
  717. sctl::StaticArray<Real, COORD_DIM*COORD_DIM> d2G = {
  718. scal*(3*dR[0]*dR[0]*invR5-invR3), scal*(3*dR[0]*dR[1]*invR5 ), scal*(3*dR[0]*dR[2]*invR5 ),
  719. scal*(3*dR[1]*dR[0]*invR5 ), scal*(3*dR[1]*dR[1]*invR5-invR3), scal*(3*dR[1]*dR[2]*invR5 ),
  720. scal*(3*dR[2]*dR[0]*invR5 ), scal*(3*dR[2]*dR[1]*invR5 ), scal*(3*dR[2]*dR[2]*invR5-invR3)
  721. };
  722. Real Un = dG[0]*Xn[0] + dG[1]*Xn[1] + dG[2]*Xn[2];
  723. Real Unn = 0;
  724. Unn += d2G[0*COORD_DIM+0]*Xn[0]*Xn[0];
  725. Unn += d2G[0*COORD_DIM+1]*Xn[0]*Xn[1];
  726. Unn += d2G[0*COORD_DIM+2]*Xn[0]*Xn[2];
  727. Unn += d2G[1*COORD_DIM+0]*Xn[1]*Xn[0];
  728. Unn += d2G[1*COORD_DIM+1]*Xn[1]*Xn[1];
  729. Unn += d2G[1*COORD_DIM+2]*Xn[1]*Xn[2];
  730. Unn += d2G[2*COORD_DIM+0]*Xn[2]*Xn[0];
  731. Unn += d2G[2*COORD_DIM+1]*Xn[2]*Xn[1];
  732. Unn += d2G[2*COORD_DIM+2]*Xn[2]*Xn[2];
  733. Real H=0;
  734. { // TODO: this is wrong
  735. Real ds0 = (dX[0*N+i]*dX[0*N+i] + dX[2*N+i]*dX[2*N+i] + dX[4*N+i]*dX[4*N+i]);
  736. Real ds1 = (dX[1*N+i]*dX[1*N+i] + dX[3*N+i]*dX[3*N+i] + dX[5*N+i]*dX[5*N+i]);
  737. H -= d2X[ 0*N+i]*Xn[0]*0.5/ds0;
  738. H -= d2X[ 4*N+i]*Xn[1]*0.5/ds0;
  739. H -= d2X[ 8*N+i]*Xn[2]*0.5/ds0;
  740. H -= d2X[ 3*N+i]*Xn[0]*0.5/ds1;
  741. H -= d2X[ 7*N+i]*Xn[1]*0.5/ds1;
  742. H -= d2X[11*N+i]*Xn[2]*0.5/ds1;
  743. }
  744. f0[i] = -2*H*Un - Unn;
  745. }
  746. }
  747. { // Compute f0 numerically
  748. sctl::Vector<Real> f1, dX1, u1;
  749. sctl::Long upsample = 3;
  750. sctl::Long Nt1 = Nt_ * upsample, Np1 = Np_ * upsample;
  751. Upsample(dX, Nt_, Np_, dX1, Nt1, Np1);
  752. Upsample(u0, Nt_, Np_, u1, Nt1, Np1);
  753. SurfaceOp<Real> SOp(comm_, Nt1, Np1);
  754. SOp.SurfLap(f1, dX1, u1);
  755. Upsample(f1, Nt1, Np1, f0, Nt_, Np_);
  756. }
  757. }
  758. }