cheb_utils.txx 32 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288
  1. /**
  2. * \file cheb_utils.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 2-11-2011
  5. * \brief This file contains chebyshev related functions.
  6. */
  7. #include <omp.h>
  8. #include <cmath>
  9. #include <cassert>
  10. #include <cstdlib>
  11. #include <iostream>
  12. #include <algorithm>
  13. #include <legendre_rule.hpp>
  14. #include <mem_utils.hpp>
  15. #include <mat_utils.hpp>
  16. #include <mem_mgr.hpp>
  17. #include <matrix.hpp>
  18. #include <profile.hpp>
  19. namespace pvfmm{
  20. template <class T>
  21. T machine_eps(){
  22. T eps=1.0;
  23. while(eps+(T)1.0>1.0) eps*=0.5;
  24. return eps;
  25. }
  26. /**
  27. * \brief Returns the values of all chebyshev polynomials up to degree d,
  28. * evaluated at points in the input vector. Output format:
  29. * { T0[in[0]], ..., T0[in[n-1]], T1[in[0]], ..., T(d-1)[in[n-1]] }
  30. */
  31. template <class T>
  32. inline void cheb_poly(int d, T* in, int n, T* out){
  33. if(d==0){
  34. for(int i=0;i<n;i++)
  35. out[i]=(fabs(in[i])<=1?1.0:0);
  36. }else if(d==1){
  37. for(int i=0;i<n;i++){
  38. out[i]=(fabs(in[i])<=1?1.0:0);
  39. out[i+n]=(fabs(in[i])<=1?in[i]:0);
  40. }
  41. }else{
  42. for(int j=0;j<n;j++){
  43. T x=(fabs(in[j])<=1?in[j]:0);
  44. T y0=(fabs(in[j])<=1?1.0:0);
  45. out[j]=y0;
  46. out[j+n]=x;
  47. T y1=x;
  48. T* y2=&out[2*n+j];
  49. for(int i=2;i<=d;i++){
  50. *y2=2*x*y1-y0;
  51. y0=y1;
  52. y1=*y2;
  53. y2=&y2[n];
  54. }
  55. }
  56. }
  57. }
  58. /**
  59. * \brief Returns the sum of the absolute value of coeffecients of the highest
  60. * order polynomial as an estimate of error.
  61. */
  62. template <class T>
  63. T cheb_err(T* cheb_coeff, int deg, int dof){
  64. T err=0;
  65. int indx=0;
  66. for(int l=0;l<dof;l++)
  67. for(int i=0;i<=deg;i++)
  68. for(int j=0;i+j<=deg;j++)
  69. for(int k=0;i+j+k<=deg;k++){
  70. if(i+j+k==deg) err+=fabs(cheb_coeff[indx]);
  71. indx++;
  72. }
  73. return err;
  74. }
  75. /**
  76. * \brief Computes Chebyshev approximation from function values at cheb node points.
  77. */
  78. template <class T, class Y>
  79. T cheb_approx(T* fn_v, int cheb_deg, int dof, T* out){
  80. //T eps=machine_eps<T>()*64;
  81. int d=cheb_deg+1;
  82. static std::vector<Matrix<Y> > precomp;
  83. static std::vector<Matrix<Y> > precomp_;
  84. Matrix<Y>* Mp ;
  85. Matrix<Y>* Mp_;
  86. #pragma omp critical (CHEB_APPROX)
  87. {
  88. if(precomp.size()<=(size_t)d){
  89. precomp .resize(d+1);
  90. precomp_.resize(d+1);
  91. }
  92. if(precomp [d].Dim(0)==0 && precomp [d].Dim(1)==0){
  93. std::vector<Y> x(d);
  94. for(int i=0;i<d;i++)
  95. x[i]=-cos((i+(T)0.5)*const_pi<T>()/d);
  96. std::vector<Y> p(d*d);
  97. cheb_poly(d-1,&x[0],d,&p[0]);
  98. for(int i=0;i<d*d;i++)
  99. p[i]=p[i]*2.0/d;
  100. Matrix<Y> Mp1(d,d,&p[0],false);
  101. Matrix<Y> Mp1_=Mp1.Transpose();
  102. precomp [d]=Mp1 ;
  103. precomp_[d]=Mp1_;
  104. }
  105. Mp =&precomp [d];
  106. Mp_=&precomp_[d];
  107. }
  108. std::vector<Y> fn_v0(d*d*d*dof);
  109. std::vector<Y> fn_v1(d*d*d);
  110. std::vector<Y> fn_v2(d*d*d);
  111. std::vector<Y> fn_v3(d*d*d);
  112. for(size_t i=0;i<(size_t)(d*d*d*dof);i++)
  113. fn_v0[i]=fn_v[i];
  114. int indx=0;
  115. for(int l=0;l<dof;l++){
  116. {
  117. Matrix<Y> M0(d*d,d,&fn_v0[d*d*d*l],false);
  118. Matrix<Y> M1(d*d,d,&fn_v1[0],false);
  119. M1=M0*(*Mp_);
  120. }
  121. {
  122. Matrix<Y> M0(d,d*d,&fn_v1[0],false);
  123. Matrix<Y> M1(d,d*d,&fn_v2[0],false);
  124. M1=(*Mp)*M0;
  125. }
  126. for(int i=0;i<d;i++){
  127. Matrix<Y> M0(d,d,&fn_v2[d*d*i],false);
  128. Matrix<Y> M1(d,d,&fn_v3[d*d*i],false);
  129. M1=(*Mp)*M0;
  130. }
  131. for(int i=0;i<d;i++)
  132. for(int j=0;j<d;j++){
  133. fn_v3[i*d+j*d*d]/=2.0;
  134. fn_v3[i+j*d*d]/=2.0;
  135. fn_v3[i+j*d]/=2.0;
  136. }
  137. Y sum=0;
  138. for(int i=0;i<d;i++)
  139. for(int j=0;i+j<d;j++)
  140. for(int k=0;i+j+k<d;k++){
  141. sum+=fabs(fn_v3[k+(j+i*d)*d]);
  142. }
  143. for(int i=0;i<d;i++)
  144. for(int j=0;i+j<d;j++)
  145. for(int k=0;i+j+k<d;k++){
  146. out[indx]=fn_v3[k+(j+i*d)*d];
  147. //if(fabs(out[indx])<eps*sum) out[indx]=0;
  148. indx++;
  149. }
  150. }
  151. return cheb_err(out,d-1,dof);
  152. }
  153. /**
  154. * \brief Returns the values of all legendre polynomials up to degree d,
  155. * evaluated at points in the input vector. Output format:
  156. * { P0[in[0]], ..., P0[in[n-1]], P1[in[0]], ..., P(d-1)[in[n-1]] }
  157. */
  158. template <class T>
  159. inline void legn_poly(int d, T* in, int n, T* out){
  160. if(d==0){
  161. for(int i=0;i<n;i++)
  162. out[i]=(fabs(in[i])<=1?1.0:0);
  163. }else if(d==1){
  164. for(int i=0;i<n;i++){
  165. out[i]=(fabs(in[i])<=1?1.0:0);
  166. out[i+n]=(fabs(in[i])<=1?in[i]:0);
  167. }
  168. }else{
  169. for(int j=0;j<n;j++){
  170. T x=(fabs(in[j])<=1?in[j]:0);
  171. T y0=(fabs(in[j])<=1?1.0:0);
  172. out[j]=y0;
  173. out[j+n]=x;
  174. T y1=x;
  175. T* y2=&out[2*n+j];
  176. for(int i=2;i<=d;i++){
  177. *y2=( (2*i-1)*x*y1-(i-1)*y0 )/i;
  178. y0=y1;
  179. y1=*y2;
  180. y2=&y2[n];
  181. }
  182. }
  183. }
  184. }
  185. /**
  186. * \brief Computes Legendre-Gauss-Lobatto nodes and weights.
  187. */
  188. template <class T>
  189. void gll_quadrature(int deg, T* x_, T* w){//*
  190. T eps=machine_eps<T>()*64;
  191. int d=deg+1;
  192. assert(d>1);
  193. int N=d-1;
  194. Vector<T> x(d,x_,false);
  195. for(int i=0;i<d;i++)
  196. x[i]=-cos((const_pi<T>()*i)/N);
  197. Matrix<T> P(d,d); P.SetZero();
  198. T err=1;
  199. Vector<T> xold(d);
  200. while(err>eps){
  201. xold=x;
  202. for(int i=0;i<d;i++){
  203. P[i][0]=1;
  204. P[i][1]=x[i];
  205. }
  206. for(int k=2;k<=N;k++)
  207. for(int i=0;i<d;i++)
  208. P[i][k]=( (2*k-1)*x[i]*P[i][k-1]-(k-1)*P[i][k-2] )/k;
  209. err=0;
  210. for(int i=0;i<d;i++){
  211. T dx=-( x[i]*P[i][N]-P[i][N-1] )/( d*P[i][N] );
  212. err=(err<fabs(dx)?fabs(dx):err);
  213. x[i]=xold[i]+dx;
  214. }
  215. }
  216. for(int i=0;i<d;i++)
  217. w[i]=2.0/(N*d*P[i][N]*P[i][N]);
  218. }
  219. /**
  220. * \brief Computes Chebyshev approximation from function values at GLL points.
  221. */
  222. template <class T, class Y>
  223. T gll2cheb(T* fn_v, int deg, int dof, T* out){//*
  224. //T eps=machine_eps<T>()*64;
  225. int d=deg+1;
  226. static std::vector<Matrix<Y> > precomp;
  227. static std::vector<Matrix<Y> > precomp_;
  228. Matrix<Y>* Mp ;
  229. Matrix<Y>* Mp_;
  230. #pragma omp critical (GLL_TO_CHEB)
  231. {
  232. if(precomp.size()<=(size_t)d){
  233. precomp .resize(d+1);
  234. precomp_.resize(d+1);
  235. std::vector<Y> x(d); //Cheb nodes.
  236. for(int i=0;i<d;i++)
  237. x[i]=-cos((i+(T)0.5)*const_pi<Y>()/d);
  238. Vector<T> w(d);
  239. Vector<T> x_legn(d); // GLL nodes.
  240. gll_quadrature(d-1, &x_legn[0], &w[0]);
  241. Matrix<T> P(d,d); //GLL node 2 GLL coeff.
  242. legn_poly(d-1,&x_legn[0],d,&P[0][0]);
  243. for(int i=0;i<d;i++)
  244. for(int j=0;j<d;j++)
  245. P[i][j]*=w[j]*0.5*(i<d-1?(2*i+1):(i));
  246. Matrix<T> M_gll2cheb(d,d); //GLL coeff 2 cheb node.
  247. legn_poly(d-1,&x[0],d,&M_gll2cheb[0][0]);
  248. Matrix<T> M_g2c; //GLL node to cheb node.
  249. M_g2c=M_gll2cheb.Transpose()*P;
  250. std::vector<Y> p(d*d);
  251. cheb_poly(d-1,&x[0],d,&p[0]);
  252. for(int i=0;i<d*d;i++)
  253. p[i]=p[i]*2.0/d;
  254. Matrix<Y> Mp1(d,d,&p[0],false);
  255. Mp1=Mp1*M_g2c;
  256. Matrix<Y> Mp1_=Mp1.Transpose();
  257. precomp [d]=Mp1 ;
  258. precomp_[d]=Mp1_;
  259. }
  260. Mp =&precomp [d];
  261. Mp_=&precomp_[d];
  262. }
  263. std::vector<Y> fn_v0(d*d*d*dof);
  264. std::vector<Y> fn_v1(d*d*d);
  265. std::vector<Y> fn_v2(d*d*d);
  266. std::vector<Y> fn_v3(d*d*d);
  267. for(size_t i=0;i<(size_t)(d*d*d*dof);i++)
  268. fn_v0[i]=fn_v[i];
  269. int indx=0;
  270. for(int l=0;l<dof;l++){
  271. {
  272. Matrix<Y> M0(d*d,d,&fn_v0[d*d*d*l],false);
  273. Matrix<Y> M1(d*d,d,&fn_v1[0],false);
  274. M1=M0*(*Mp_);
  275. }
  276. {
  277. Matrix<Y> M0(d,d*d,&fn_v1[0],false);
  278. Matrix<Y> M1(d,d*d,&fn_v2[0],false);
  279. M1=(*Mp)*M0;
  280. }
  281. for(int i=0;i<d;i++){
  282. Matrix<Y> M0(d,d,&fn_v2[d*d*i],false);
  283. Matrix<Y> M1(d,d,&fn_v3[d*d*i],false);
  284. M1=(*Mp)*M0;
  285. }
  286. for(int i=0;i<d;i++)
  287. for(int j=0;j<d;j++){
  288. fn_v3[i*d+j*d*d]/=2.0;
  289. fn_v3[i+j*d*d]/=2.0;
  290. fn_v3[i+j*d]/=2.0;
  291. }
  292. Y sum=0;
  293. for(int i=0;i<d;i++)
  294. for(int j=0;i+j<d;j++)
  295. for(int k=0;i+j+k<d;k++){
  296. sum+=fabs(fn_v3[k+(j+i*d)*d]);
  297. }
  298. for(int i=0;i<d;i++)
  299. for(int j=0;i+j<d;j++)
  300. for(int k=0;i+j+k<d;k++){
  301. out[indx]=fn_v3[k+(j+i*d)*d];
  302. //if(fabs(out[indx])<eps*sum) out[indx]=0;
  303. indx++;
  304. }
  305. }
  306. return cheb_err(out,d-1,dof);
  307. }
  308. /**
  309. * \brief Computes Chebyshev approximation from the input function pointer.
  310. */
  311. template <class T>
  312. T cheb_approx(T (*fn)(T,T,T), int cheb_deg, T* coord, T s, std::vector<T>& out){
  313. int d=cheb_deg+1;
  314. std::vector<T> x(d);
  315. for(int i=0;i<d;i++)
  316. x[i]=cos((i+(T)0.5)*const_pi<T>()/d);
  317. std::vector<T> p;
  318. cheb_poly(d-1,&x[0],d,&p[0]);
  319. std::vector<T> x1(d);
  320. std::vector<T> x2(d);
  321. std::vector<T> x3(d);
  322. for(int i=0;i<d;i++){
  323. x1[i]=(x[i]+1.0)/2.0*s+coord[0];
  324. x2[i]=(x[i]+1.0)/2.0*s+coord[1];
  325. x3[i]=(x[i]+1.0)/2.0*s+coord[2];
  326. }
  327. std::vector<T> fn_v(d*d*d);
  328. T* fn_p=&fn_v[0];
  329. for(int i=0;i<d;i++){
  330. for(int j=0;j<d;j++){
  331. for(int k=0;k<d;k++){
  332. *fn_p=fn(x3[k],x2[j],x1[i]);
  333. fn_p++;
  334. }
  335. }
  336. }
  337. out.resize((d*(d+1)*(d+2))/6);
  338. return cheb_approx(&fn_v[0], d-1, 1, &out[0]);
  339. }
  340. /**
  341. * \brief Evaluates polynomial values from input coefficients at points on
  342. * a regular grid defined by in_x, in_y, in_z the values in the input vector.
  343. */
  344. template <class T>
  345. void cheb_eval(Vector<T>& coeff_, int cheb_deg, std::vector<T>& in_x, std::vector<T>& in_y, std::vector<T>& in_z, Vector<T>& out){
  346. int d=cheb_deg+1;
  347. int dof=coeff_.Dim()/((d*(d+1)*(d+2))/6);
  348. assert(coeff_.Dim()==(size_t)(d*(d+1)*(d+2)*dof)/6);
  349. std::vector<T> coeff(d*d*d*dof);
  350. {// Rearrange data
  351. int indx=0;
  352. for(int l=0;l<dof;l++)
  353. for(int i=0;i<d;i++)
  354. for(int j=0;i+j<d;j++)
  355. for(int k=0;i+j+k<d;k++){
  356. coeff[(k+(j+(i+l*d)*d)*d)]=coeff_[indx];
  357. indx++;
  358. }
  359. }
  360. int n1=in_x.size();
  361. int n2=in_y.size();
  362. int n3=in_z.size();
  363. out.Resize(n1*n2*n3*dof);
  364. if(n1==0 || n2==0 || n3==0) return;
  365. std::vector<T> p1(n1*d);
  366. std::vector<T> p2(n2*d);
  367. std::vector<T> p3(n3*d);
  368. cheb_poly(d-1,&in_x[0],n1,&p1[0]);
  369. cheb_poly(d-1,&in_y[0],n2,&p2[0]);
  370. cheb_poly(d-1,&in_z[0],n3,&p3[0]);
  371. std::vector<T> fn_v1(n1*d *d );
  372. std::vector<T> fn_v2(n1*d *n3);
  373. std::vector<T> fn_v3(n1*n2*n3);
  374. Matrix<T> Mp1(d,n1,&p1[0],false);
  375. Matrix<T> Mp2(d,n2,&p2[0],false);
  376. Matrix<T> Mp3(d,n3,&p3[0],false);
  377. Matrix<T> Mp2_=Mp2.Transpose();
  378. Matrix<T> Mp3_=Mp3.Transpose();
  379. for(int k=0;k<dof;k++){
  380. {
  381. Matrix<T> M0(d*d,d,&coeff[k*d*d*d],false);
  382. Matrix<T> M1(d*d,n1,&fn_v1[0],false);
  383. M1=M0*Mp1;
  384. }
  385. {
  386. Matrix<T> M0(d,d*n1,&fn_v1[0],false);
  387. Matrix<T> M1(n3,d*n1,&fn_v2[0],false);
  388. M1=Mp3_*M0;
  389. }
  390. {
  391. int dn1=d*n1;
  392. int n2n1=n2*n1;
  393. for(int i=0;i<n3;i++){
  394. Matrix<T> M0(d,n1,&fn_v2[i*dn1],false);
  395. Matrix<T> M1(n2,n1,&fn_v3[i*n2n1],false);
  396. M1=Mp2_*M0;
  397. }
  398. }
  399. mem::memcopy(&out[n1*n2*n3*k],&fn_v3[0],n1*n2*n3*sizeof(T));
  400. }
  401. }
  402. /**
  403. * \brief Evaluates polynomial values from input coefficients at points
  404. * in the coord vector.
  405. */
  406. template <class T>
  407. inline void cheb_eval(Vector<T>& coeff_, int cheb_deg, std::vector<T>& coord, Vector<T>& out){
  408. int dim=3;
  409. int d=cheb_deg+1;
  410. int n=coord.size()/dim;
  411. int dof=coeff_.Dim()/((d*(d+1)*(d+2))/6);
  412. assert(coeff_.Dim()==(size_t)(d*(d+1)*(d+2)*dof)/6);
  413. std::vector<T> coeff(d*d*d*dof);
  414. {// Rearrange data
  415. int indx=0;
  416. for(int l=0;l<dof;l++)
  417. for(int i=0;i<d;i++)
  418. for(int j=0;i+j<d;j++)
  419. for(int k=0;i+j+k<d;k++){
  420. coeff[(k+(j+(i+l*d)*d)*d)]=coeff_[indx];
  421. indx++;
  422. }
  423. }
  424. Matrix<T> coord_(n,dim,&coord[0]);
  425. coord_=coord_.Transpose();
  426. Matrix<T> px(d,n);
  427. Matrix<T> py(d,n);
  428. Matrix<T> pz(d,n);
  429. cheb_poly(d-1,&(coord_[0][0]),n,&(px[0][0]));
  430. cheb_poly(d-1,&(coord_[1][0]),n,&(py[0][0]));
  431. cheb_poly(d-1,&(coord_[2][0]),n,&(pz[0][0]));
  432. Matrix<T> M_coeff0(d*d*dof, d, &coeff[0], false);
  433. Matrix<T> M0 = (M_coeff0 * px).Transpose(); // {n, dof*d*d}
  434. py = py.Transpose();
  435. pz = pz.Transpose();
  436. out.Resize(n*dof);
  437. for(int i=0; i<n; i++)
  438. for(int j=0; j<dof; j++){
  439. Matrix<T> M0_ (d, d, &(M0[i][ j*d*d]), false);
  440. Matrix<T> py_ (d, 1, &(py[i][ 0]), false);
  441. Matrix<T> pz_ (1, d, &(pz[i][ 0]), false);
  442. Matrix<T> M_out(1, 1, &( out[i*dof+j]), false);
  443. M_out += pz_ * M0_ * py_;
  444. }
  445. }
  446. /**
  447. * \brief Returns the values of all Chebyshev basis functions of degree up to d
  448. * evaluated at the point coord.
  449. */
  450. template <class T>
  451. inline void cheb_eval(int cheb_deg, T* coord, T* coeff0,T* buff){
  452. int d=cheb_deg+1;
  453. std::vector<T> coeff(d*d*d);
  454. T* p=&buff[0];
  455. T* p_=&buff[3*d];
  456. cheb_poly(d-1,&coord[0],3,&p[0]);
  457. for(int i=0;i<d;i++){
  458. p_[i]=p[i*3];
  459. p_[i+d]=p[i*3+1];
  460. p_[i+2*d]=p[i*3+2];
  461. }
  462. T* coeff_=&buff[2*3*d];
  463. Matrix<T> v_p0 (1, d, & p_[0],false);
  464. Matrix<T> v_p1 (d, 1, & p_[d],false);
  465. Matrix<T> M_coeff_(d, d, &coeff_[0],false);
  466. M_coeff_ = v_p1 * v_p0; // */
  467. //mat::gemm(CblasRowMajor,CblasNoTrans,CblasNoTrans,d,d,1,1.0,&p_[d],1,&p_[0],d,0.0,&coeff_[0],d);
  468. Matrix<T> v_p2 (d, 1, & p_[2*d],false);
  469. Matrix<T> v_coeff_(1, d*d, &coeff_[ 0],false);
  470. Matrix<T> M_coeff (d, d*d, &coeff [ 0],false);
  471. M_coeff = v_p2 * v_coeff_; // */
  472. //mat::gemm(CblasRowMajor,CblasNoTrans,CblasNoTrans,d,d*d,1,1.0,&p_[2*d],1,&coeff_[0],d*d,0.0,&coeff[0],d*d);
  473. {// Rearrange data
  474. int indx=0;
  475. for(int i=0;i<d;i++)
  476. for(int j=0;i+j<d;j++)
  477. for(int k=0;i+j+k<d;k++){
  478. coeff0[indx]=coeff[(k+(j+i*d)*d)];
  479. indx++;
  480. }
  481. }
  482. }
  483. /**
  484. * \brief Computes a least squares solution for Chebyshev approximation over a
  485. * cube from point samples.
  486. * \param[in] deg Maximum degree of the polynomial.
  487. * \param[in] coord Coordinates of points (x,y,z interleaved).
  488. * \param[in] node_coord Coordinates of the octant.
  489. * \param[in] node_size Length of the side of the octant.
  490. * \param[out] cheb_coeff Output coefficients.
  491. */
  492. template <class T>
  493. void points2cheb(int deg, T* coord, T* val, int n, int dim, T* node_coord, T node_size, Vector<T>& cheb_coeff){
  494. if(n==0) return;
  495. int deg_=((int)(pow((T)n*6,1.0/3.0)+0.5))/2;
  496. deg_=(deg_>deg?deg:deg_);
  497. deg_=(deg_>0?deg_:1);
  498. int deg3=((deg_+1)*(deg_+2)*(deg_+3))/6;
  499. cheb_coeff.Resize(dim*((deg+1)*(deg+2)*(deg+3))/6);
  500. cheb_coeff.SetZero();
  501. //Map coordinates to unit cube
  502. std::vector<T> coord_(n*3);
  503. for(int i=0;i<n;i++){
  504. coord_[i*3 ]=(coord[i*3 ]-node_coord[0])*2.0/node_size-1.0;
  505. coord_[i*3+1]=(coord[i*3+1]-node_coord[1])*2.0/node_size-1.0;
  506. coord_[i*3+2]=(coord[i*3+2]-node_coord[2])*2.0/node_size-1.0;
  507. }
  508. //Compute the matrix M
  509. Matrix<T> M(n,deg3);
  510. std::vector<T> buff((deg_+1)*(deg_+1+3*2));
  511. for(int i=0;i<n;i++)
  512. cheb_eval(deg_,&coord_[i*3],&(M[i][0]),&buff[0]);
  513. //Compute the pinv and get the cheb_coeff.
  514. Matrix<T> M_val(n,dim,&val[0]);
  515. T eps=machine_eps<T>()*64;
  516. Matrix<T> cheb_coeff_=(M.pinv(eps)*M_val).Transpose();
  517. //Set the output
  518. int indx=0;
  519. int indx1=0;
  520. for(int l=0;l<dim;l++)
  521. for(int i=0;i <=deg;i++)
  522. for(int j=0;i+j <=deg;j++)
  523. for(int k=0;i+j+k<=deg;k++){
  524. if(i+j+k<=deg_){
  525. cheb_coeff[indx]=cheb_coeff_[0][indx1];
  526. indx1++;
  527. }else{
  528. cheb_coeff[indx]=0;
  529. }
  530. indx++;
  531. }
  532. }
  533. template <class T>
  534. void quad_rule(int n, T* x, T* w){//*
  535. static std::vector<Vector<T> > x_lst(10000);
  536. static std::vector<Vector<T> > w_lst(10000);
  537. assert(n<10000);
  538. bool done=false;
  539. #pragma omp critical (QUAD_RULE)
  540. if(x_lst[n].Dim()>0){
  541. Vector<T>& x_=x_lst[n];
  542. Vector<T>& w_=w_lst[n];
  543. for(int i=0;i<n;i++){
  544. x[i]=x_[i];
  545. w[i]=w_[i];
  546. }
  547. done=true;
  548. }
  549. if(done) return;
  550. Vector<T> x_(n);
  551. Vector<T> w_(n);
  552. { //Gauss-Chebyshev quadrature nodes and weights
  553. for(int i=0;i<n;i++){
  554. x_[i]=-cos((T)(2.0*i+1.0)/(2.0*n)*const_pi<T>());
  555. w_[i]=0;//sqrt(1.0-x_[i]*x_[i])*const_pi<T>()/n;
  556. }
  557. Matrix<T> M(n,n);
  558. cheb_poly(n-1, &x_[0], n, &M[0][0]);
  559. for(size_t i=0;i<n;i++) M[0][i]/=2.0;
  560. std::vector<T> w_sample(n,0);
  561. if(n>0) w_sample[0]=2.0;
  562. if(n>1) w_sample[1]=0.0;
  563. if(n>2) w_sample[2]=-((T)2.0)/3;
  564. if(n>3) w_sample[3]=0.0;
  565. if(n>4) w_sample[4]=-((T)2.0)/15;
  566. if(n>5) w_sample[5]=0.0;
  567. if(n>6) w_sample[5]=((T)64)/7-((T)96)/5+((T)36)/3-2;
  568. if(n>7) w_sample[5]=0;
  569. if(n>8){
  570. T eps=machine_eps<T>()*64;
  571. std::vector<T> qx(n-1);
  572. std::vector<T> qw(n-1);
  573. quad_rule(n-1, &qx[0], &qw[0]);
  574. T err=1.0;
  575. std::vector<T> w_prev;
  576. for(size_t iter=1;err>eps*iter;iter*=2){
  577. w_prev=w_sample;
  578. w_sample.assign(n,0);
  579. size_t N=(n-1)*iter;
  580. std::vector<T> x_sample(N,0);
  581. Matrix<T> M_sample(n,N);
  582. for(size_t i=0;i<iter;i++){
  583. for(size_t j=0;j<n-1;j++){
  584. x_sample[j+i*(n-1)]=(2*i+qx[j]+1)/iter-1;
  585. }
  586. }
  587. cheb_poly(n-1, &x_sample[0], N, &M_sample[0][0]);
  588. for(size_t i=0;i<n;i++)
  589. for(size_t j=0;j<iter;j++)
  590. for(size_t k=0;k<n-1;k++){
  591. w_sample[i]+=M_sample[i][k+j*(n-1)]*qw[k];
  592. }
  593. for(size_t i=0;i<n;i++) w_sample[i]/=iter;
  594. for(size_t i=1;i<n;i+=2) w_sample[i]=0.0;
  595. err=0;
  596. for(size_t i=0;i<n;i++) err+=fabs(w_sample[i]-w_prev[i]);
  597. }
  598. }
  599. for(size_t i=0;i<n;i++)
  600. for(size_t j=0;j<n;j++){
  601. M[i][j]*=w_sample[i];
  602. }
  603. for(size_t i=0;i<n;i++)
  604. for(size_t j=0;j<n;j++){
  605. w_[j]+=M[i][j]*2/n;
  606. }
  607. }
  608. { //Trapezoidal quadrature nodes and weights
  609. //for(int i=0;i<n;i++){
  610. // x_[i]=(2.0*i+1.0)/(1.0*n)-1.0;
  611. // w_[i]=2.0/n;
  612. //}
  613. }
  614. #pragma omp critical (QUAD_RULE)
  615. { // Set x_lst, w_lst
  616. x_lst[n]=x_;
  617. w_lst[n]=w_;
  618. }
  619. quad_rule(n, x, w);
  620. }
  621. template <>
  622. void quad_rule<double>(int n, double* x, double* w){//*
  623. static std::vector<Vector<double> > x_lst(10000);
  624. static std::vector<Vector<double> > w_lst(10000);
  625. assert(n<10000);
  626. bool done=false;
  627. #pragma omp critical (QUAD_RULE)
  628. if(x_lst[n].Dim()>0){
  629. Vector<double>& x_=x_lst[n];
  630. Vector<double>& w_=w_lst[n];
  631. for(int i=0;i<n;i++){
  632. x[i]=x_[i];
  633. w[i]=w_[i];
  634. }
  635. done=true;
  636. }
  637. if(done) return;
  638. Vector<double> x_(n);
  639. Vector<double> w_(n);
  640. { //Gauss-Legendre quadrature nodes and weights
  641. double alpha=0.0;
  642. double beta=0.0;
  643. double a=-1.0;
  644. double b= 1.0;
  645. int kind = 1;
  646. cgqf ( n, kind, (double)alpha, (double)beta, (double)a, (double)b, &x_[0], &w_[0] );
  647. }
  648. #pragma omp critical (QUAD_RULE)
  649. { // Set x_lst, w_lst
  650. x_lst[n]=x_;
  651. w_lst[n]=w_;
  652. }
  653. quad_rule(n, x, w);
  654. }
  655. template <class T>
  656. std::vector<T> integ_pyramid(int m, T* s, T r, int nx, Kernel<T>& kernel, int* perm){//*
  657. static mem::MemoryManager mem_mgr(16*1024*1024*sizeof(T));
  658. int ny=nx;
  659. int nz=nx;
  660. T eps=machine_eps<T>()*64;
  661. int k_dim=kernel.ker_dim[0]*kernel.ker_dim[1];
  662. std::vector<T> qp_x(nx), qw_x(nx);
  663. std::vector<T> qp_y(ny), qw_y(ny);
  664. std::vector<T> qp_z(nz), qw_z(nz);
  665. std::vector<T> p_x(nx*m);
  666. std::vector<T> p_y(ny*m);
  667. std::vector<T> p_z(nz*m);
  668. std::vector<T> x_;
  669. { // Build stack along X-axis.
  670. x_.push_back(s[0]);
  671. x_.push_back(fabs(1.0-s[0])+s[0]);
  672. x_.push_back(fabs(1.0-s[1])+s[0]);
  673. x_.push_back(fabs(1.0+s[1])+s[0]);
  674. x_.push_back(fabs(1.0-s[2])+s[0]);
  675. x_.push_back(fabs(1.0+s[2])+s[0]);
  676. std::sort(x_.begin(),x_.end());
  677. for(int i=0;i<x_.size();i++){
  678. if(x_[i]<-1.0) x_[i]=-1.0;
  679. if(x_[i]> 1.0) x_[i]= 1.0;
  680. }
  681. std::vector<T> x_new;
  682. T x_jump=fabs(1.0-s[0]);
  683. if(fabs(1.0-s[1])>eps) x_jump=std::min(x_jump,(T)fabs(1.0-s[1]));
  684. if(fabs(1.0+s[1])>eps) x_jump=std::min(x_jump,(T)fabs(1.0+s[1]));
  685. if(fabs(1.0-s[2])>eps) x_jump=std::min(x_jump,(T)fabs(1.0-s[2]));
  686. if(fabs(1.0+s[2])>eps) x_jump=std::min(x_jump,(T)fabs(1.0+s[2]));
  687. for(int k=0; k<x_.size()-1; k++){
  688. T x0=x_[k];
  689. T x1=x_[k+1];
  690. T A0=0;
  691. T A1=0;
  692. { // A0
  693. T y0=s[1]-(x0-s[0]); if(y0<-1.0) y0=-1.0; if(y0> 1.0) y0= 1.0;
  694. T y1=s[1]+(x0-s[0]); if(y1<-1.0) y1=-1.0; if(y1> 1.0) y1= 1.0;
  695. T z0=s[2]-(x0-s[0]); if(z0<-1.0) z0=-1.0; if(z0> 1.0) z0= 1.0;
  696. T z1=s[2]+(x0-s[0]); if(z1<-1.0) z1=-1.0; if(z1> 1.0) z1= 1.0;
  697. A0=(y1-y0)*(z1-z0);
  698. }
  699. { // A1
  700. T y0=s[1]-(x1-s[0]); if(y0<-1.0) y0=-1.0; if(y0> 1.0) y0= 1.0;
  701. T y1=s[1]+(x1-s[0]); if(y1<-1.0) y1=-1.0; if(y1> 1.0) y1= 1.0;
  702. T z0=s[2]-(x1-s[0]); if(z0<-1.0) z0=-1.0; if(z0> 1.0) z0= 1.0;
  703. T z1=s[2]+(x1-s[0]); if(z1<-1.0) z1=-1.0; if(z1> 1.0) z1= 1.0;
  704. A1=(y1-y0)*(z1-z0);
  705. }
  706. T V=0.5*(A0+A1)*(x1-x0);
  707. if(V<eps) continue;
  708. if(!x_new.size()) x_new.push_back(x0);
  709. x_jump=std::max(x_jump,x0-s[0]);
  710. while(s[0]+x_jump*1.5<x1){
  711. x_new.push_back(s[0]+x_jump);
  712. x_jump*=2.0;
  713. }
  714. if(x_new.back()+eps<x1) x_new.push_back(x1);
  715. }
  716. assert(x_new.size()<30);
  717. x_.swap(x_new);
  718. }
  719. Vector<T> k_out( ny*nz*k_dim,(T*)mem_mgr.malloc( ny*nz*k_dim*sizeof(T)),false); //Output of kernel evaluation.
  720. Vector<T> I0 ( ny*m *k_dim,(T*)mem_mgr.malloc( ny*m *k_dim*sizeof(T)),false);
  721. Vector<T> I1 ( m *m *k_dim,(T*)mem_mgr.malloc( m *m *k_dim*sizeof(T)),false);
  722. Vector<T> I2 (m *m *m *k_dim,(T*)mem_mgr.malloc(m *m *m *k_dim*sizeof(T)),false); I2.SetZero();
  723. if(x_.size()>1)
  724. for(int k=0; k<x_.size()-1; k++){
  725. T x0=x_[k];
  726. T x1=x_[k+1];
  727. { // Set qp_x
  728. std::vector<T> qp(nx);
  729. std::vector<T> qw(nx);
  730. quad_rule(nx,&qp[0],&qw[0]);
  731. for(int i=0; i<nx; i++)
  732. qp_x[i]=(x1-x0)*qp[i]/2.0+(x1+x0)/2.0;
  733. qw_x=qw;
  734. }
  735. cheb_poly(m-1,&qp_x[0],nx,&p_x[0]);
  736. for(int i=0; i<nx; i++){
  737. T y0=s[1]-(qp_x[i]-s[0]); if(y0<-1.0) y0=-1.0; if(y0> 1.0) y0= 1.0;
  738. T y1=s[1]+(qp_x[i]-s[0]); if(y1<-1.0) y1=-1.0; if(y1> 1.0) y1= 1.0;
  739. T z0=s[2]-(qp_x[i]-s[0]); if(z0<-1.0) z0=-1.0; if(z0> 1.0) z0= 1.0;
  740. T z1=s[2]+(qp_x[i]-s[0]); if(z1<-1.0) z1=-1.0; if(z1> 1.0) z1= 1.0;
  741. { // Set qp_y
  742. std::vector<T> qp(ny);
  743. std::vector<T> qw(ny);
  744. quad_rule(ny,&qp[0],&qw[0]);
  745. for(int j=0; j<ny; j++)
  746. qp_y[j]=(y1-y0)*qp[j]/2.0+(y1+y0)/2.0;
  747. qw_y=qw;
  748. }
  749. { // Set qp_z
  750. std::vector<T> qp(nz);
  751. std::vector<T> qw(nz);
  752. quad_rule(nz,&qp[0],&qw[0]);
  753. for(int j=0; j<nz; j++)
  754. qp_z[j]=(z1-z0)*qp[j]/2.0+(z1+z0)/2.0;
  755. qw_z=qw;
  756. }
  757. cheb_poly(m-1,&qp_y[0],ny,&p_y[0]);
  758. cheb_poly(m-1,&qp_z[0],nz,&p_z[0]);
  759. { // k_out = kernel x qw
  760. T src[3]={0,0,0};
  761. std::vector<T> trg(ny*nz*3);
  762. for(int i0=0; i0<ny; i0++){
  763. size_t indx0=i0*nz*3;
  764. for(int i1=0; i1<nz; i1++){
  765. size_t indx1=indx0+i1*3;
  766. trg[indx1+perm[0]]=(s[0]-qp_x[i ])*r*0.5*perm[1];
  767. trg[indx1+perm[2]]=(s[1]-qp_y[i0])*r*0.5*perm[3];
  768. trg[indx1+perm[4]]=(s[2]-qp_z[i1])*r*0.5*perm[5];
  769. }
  770. }
  771. {
  772. Matrix<T> k_val(ny*nz*kernel.ker_dim[0],kernel.ker_dim[1]);
  773. kernel.BuildMatrix(&src[0],1,&trg[0],ny*nz,&k_val[0][0]);
  774. Matrix<T> k_val_t(kernel.ker_dim[1],ny*nz*kernel.ker_dim[0],&k_out[0], false);
  775. k_val_t=k_val.Transpose();
  776. }
  777. for(int kk=0; kk<k_dim; kk++){
  778. for(int i0=0; i0<ny; i0++){
  779. size_t indx=(kk*ny+i0)*nz;
  780. for(int i1=0; i1<nz; i1++){
  781. k_out[indx+i1] *= qw_y[i0]*qw_z[i1];
  782. }
  783. }
  784. }
  785. }
  786. I0.SetZero();
  787. for(int kk=0; kk<k_dim; kk++){
  788. for(int i0=0; i0<ny; i0++){
  789. size_t indx0=(kk*ny+i0)*nz;
  790. size_t indx1=(kk*ny+i0)* m;
  791. for(int i2=0; i2<m; i2++){
  792. for(int i1=0; i1<nz; i1++){
  793. I0[indx1+i2] += k_out[indx0+i1]*p_z[i2*nz+i1];
  794. }
  795. }
  796. }
  797. }
  798. I1.SetZero();
  799. for(int kk=0; kk<k_dim; kk++){
  800. for(int i2=0; i2<ny; i2++){
  801. size_t indx0=(kk*ny+i2)*m;
  802. for(int i0=0; i0<m; i0++){
  803. size_t indx1=(kk* m+i0)*m;
  804. T py=p_y[i0*ny+i2];
  805. for(int i1=0; i0+i1<m; i1++){
  806. I1[indx1+i1] += I0[indx0+i1]*py;
  807. }
  808. }
  809. }
  810. }
  811. T v=(x1-x0)*(y1-y0)*(z1-z0);
  812. for(int kk=0; kk<k_dim; kk++){
  813. for(int i0=0; i0<m; i0++){
  814. T px=p_x[i+i0*nx]*qw_x[i]*v;
  815. for(int i1=0; i0+i1<m; i1++){
  816. size_t indx0= (kk*m+i1)*m;
  817. size_t indx1=((kk*m+i0)*m+i1)*m;
  818. for(int i2=0; i0+i1+i2<m; i2++){
  819. I2[indx1+i2] += I1[indx0+i2]*px;
  820. }
  821. }
  822. }
  823. }
  824. }
  825. }
  826. for(int i=0;i<m*m*m*k_dim;i++)
  827. I2[i]=I2[i]*r*r*r/64.0;
  828. if(x_.size()>1)
  829. Profile::Add_FLOP(( 2*ny*nz*m*k_dim
  830. +ny*m*(m+1)*k_dim
  831. +2*m*(m+1)*k_dim
  832. +m*(m+1)*(m+2)/3*k_dim)*nx*(x_.size()-1));
  833. std::vector<T> I2_(&I2[0], &I2[0]+I2.Dim());
  834. mem_mgr.free(&k_out[0]);
  835. mem_mgr.free(&I0 [0]);
  836. mem_mgr.free(&I1 [0]);
  837. mem_mgr.free(&I2 [0]);
  838. return I2_;
  839. }
  840. template <class T>
  841. std::vector<T> integ(int m, T* s, T r, int n, Kernel<T>& kernel){//*
  842. //Compute integrals over pyramids in all directions.
  843. int k_dim=kernel.ker_dim[0]*kernel.ker_dim[1];
  844. T s_[3];
  845. s_[0]=s[0]*2.0/r-1.0;
  846. s_[1]=s[1]*2.0/r-1.0;
  847. s_[2]=s[2]*2.0/r-1.0;
  848. T s1[3];
  849. int perm[6];
  850. std::vector<T> U_[6];
  851. s1[0]= s_[0];s1[1]=s_[1];s1[2]=s_[2];
  852. perm[0]= 0; perm[2]= 1; perm[4]= 2;
  853. perm[1]= 1; perm[3]= 1; perm[5]= 1;
  854. U_[0]=integ_pyramid<T>(m,s1,r,n,kernel,perm);
  855. s1[0]=-s_[0];s1[1]=s_[1];s1[2]=s_[2];
  856. perm[0]= 0; perm[2]= 1; perm[4]= 2;
  857. perm[1]=-1; perm[3]= 1; perm[5]= 1;
  858. U_[1]=integ_pyramid<T>(m,s1,r,n,kernel,perm);
  859. s1[0]= s_[1];s1[1]=s_[0];s1[2]=s_[2];
  860. perm[0]= 1; perm[2]= 0; perm[4]= 2;
  861. perm[1]= 1; perm[3]= 1; perm[5]= 1;
  862. U_[2]=integ_pyramid<T>(m,s1,r,n,kernel,perm);
  863. s1[0]=-s_[1];s1[1]=s_[0];s1[2]=s_[2];
  864. perm[0]= 1; perm[2]= 0; perm[4]= 2;
  865. perm[1]=-1; perm[3]= 1; perm[5]= 1;
  866. U_[3]=integ_pyramid<T>(m,s1,r,n,kernel,perm);
  867. s1[0]= s_[2];s1[1]=s_[0];s1[2]=s_[1];
  868. perm[0]= 2; perm[2]= 0; perm[4]= 1;
  869. perm[1]= 1; perm[3]= 1; perm[5]= 1;
  870. U_[4]=integ_pyramid<T>(m,s1,r,n,kernel,perm);
  871. s1[0]=-s_[2];s1[1]=s_[0];s1[2]=s_[1];
  872. perm[0]= 2; perm[2]= 0; perm[4]= 1;
  873. perm[1]=-1; perm[3]= 1; perm[5]= 1;
  874. U_[5]=integ_pyramid<T>(m,s1,r,n,kernel,perm);
  875. std::vector<T> U; U.assign(m*m*m*k_dim,0);
  876. for(int kk=0; kk<k_dim; kk++){
  877. for(int i=0;i<m;i++){
  878. for(int j=0;j<m;j++){
  879. for(int k=0;k<m;k++){
  880. U[kk*m*m*m + k*m*m + j*m + i]+=U_[0][kk*m*m*m + i*m*m + j*m + k];
  881. U[kk*m*m*m + k*m*m + j*m + i]+=U_[1][kk*m*m*m + i*m*m + j*m + k]*(i%2?-1.0:1.0);
  882. }
  883. }
  884. }
  885. }
  886. for(int kk=0; kk<k_dim; kk++){
  887. for(int i=0; i<m; i++){
  888. for(int j=0; j<m; j++){
  889. for(int k=0; k<m; k++){
  890. U[kk*m*m*m + k*m*m + i*m + j]+=U_[2][kk*m*m*m + i*m*m + j*m + k];
  891. U[kk*m*m*m + k*m*m + i*m + j]+=U_[3][kk*m*m*m + i*m*m + j*m + k]*(i%2?-1.0:1.0);
  892. }
  893. }
  894. }
  895. }
  896. for(int kk=0; kk<k_dim; kk++){
  897. for(int i=0; i<m; i++){
  898. for(int j=0; j<m; j++){
  899. for(int k=0; k<m; k++){
  900. U[kk*m*m*m + i*m*m + k*m + j]+=U_[4][kk*m*m*m + i*m*m + j*m + k];
  901. U[kk*m*m*m + i*m*m + k*m + j]+=U_[5][kk*m*m*m + i*m*m + j*m + k]*(i%2?-1.0:1.0);
  902. }
  903. }
  904. }
  905. }
  906. return U;
  907. }
  908. /**
  909. * \brief
  910. * \param[in] r Length of the side of cubic region.
  911. */
  912. template <class T>
  913. std::vector<T> cheb_integ(int m, T* s_, T r_, Kernel<T>& kernel){
  914. T eps=machine_eps<T>();
  915. T r=r_;
  916. T s[3]={s_[0],s_[1],s_[2]};
  917. int n=m+1;
  918. T err=1.0;
  919. int k_dim=kernel.ker_dim[0]*kernel.ker_dim[1];
  920. std::vector<T> U=integ<T>(m+1,s,r,n,kernel);
  921. std::vector<T> U_;
  922. while(err>eps*n*n){
  923. n=(int)round(n*1.3);
  924. if(n>300){
  925. std::cout<<"Cheb_Integ::Failed to converge.[";
  926. std::cout<<((double)err )<<",";
  927. std::cout<<((double)s[0])<<",";
  928. std::cout<<((double)s[1])<<",";
  929. std::cout<<((double)s[2])<<"]\n";
  930. break;
  931. }
  932. U_=integ<T>(m+1,s,r,n,kernel);
  933. err=0;
  934. for(int i=0;i<(m+1)*(m+1)*(m+1)*k_dim;i++)
  935. if(fabs(U[i]-U_[i])>err)
  936. err=fabs(U[i]-U_[i]);
  937. U=U_;
  938. }
  939. std::vector<T> U0(((m+1)*(m+2)*(m+3)*k_dim)/6);
  940. {// Rearrange data
  941. int indx=0;
  942. int* ker_dim=kernel.ker_dim;
  943. for(int l0=0;l0<ker_dim[0];l0++)
  944. for(int l1=0;l1<ker_dim[1];l1++)
  945. for(int i=0;i<=m;i++)
  946. for(int j=0;i+j<=m;j++)
  947. for(int k=0;i+j+k<=m;k++){
  948. U0[indx]=U[(k+(j+(i+(l1*ker_dim[0]+l0)*(m+1))*(m+1))*(m+1))];
  949. indx++;
  950. }
  951. }
  952. return U0;
  953. }
  954. template <class T>
  955. std::vector<T> cheb_nodes(int deg, int dim){
  956. int d=deg+1;
  957. std::vector<T> x(d);
  958. for(int i=0;i<d;i++)
  959. x[i]=-cos((i+(T)0.5)*const_pi<T>()/d)*0.5+0.5;
  960. if(dim==1) return x;
  961. int n1=(int)(pow((T)d,dim)+0.5);
  962. std::vector<T> y(n1*dim);
  963. for(int i=0;i<dim;i++){
  964. int n2=(int)(pow((T)d,i)+0.5);
  965. for(int j=0;j<n1;j++){
  966. y[j*dim+i]=x[(j/n2)%d];
  967. }
  968. }
  969. return y;
  970. }
  971. template <class T>
  972. void cheb_diff(T* A, int deg, T* B){
  973. int d=deg+1;
  974. static Matrix<T> M;
  975. #pragma omp critical (CHEB_DIFF)
  976. if(M.Dim(0)!=d){
  977. M.Resize(d,d);
  978. for(int i=0;i<d;i++){
  979. for(int j=0;j<d;j++) M[j][i]=0;
  980. for(int j=1-(i%2);j<i-1;j=j+2){
  981. M[j][i]=2*i*2;
  982. }
  983. if(i%2==1) M[0][i]-=i*2;
  984. }
  985. }
  986. Matrix<T> MA(d,1,A,false);
  987. Matrix<T> MB(d,1,B,false);
  988. MB=M*MA;
  989. }
  990. template <class T>
  991. void cheb_diff(T* A, int deg, int dim, int curr_dim, T* B){
  992. int d=deg+1;
  993. static Matrix<T> M;
  994. #pragma omp critical (CHEB_DIFF1)
  995. if(M.Dim(0)!=(size_t)d){
  996. M.Resize(d,d);
  997. for(int i=0;i<d;i++){
  998. for(int j=0;j<d;j++) M[j][i]=0;
  999. for(int j=1-(i%2);j<i;j=j+2){
  1000. M[j][i]=2*i*2;
  1001. }
  1002. if(i%2==1) M[0][i]-=i*2;
  1003. }
  1004. }
  1005. int n1=(int)(pow((T)d,curr_dim)+0.5);
  1006. int n2=(int)(pow((T)d,dim-curr_dim-1)+0.5);
  1007. for(int i=0;i<n2;i++){
  1008. Matrix<T> MA(d,n1,&A[i*n1*d],false);
  1009. Matrix<T> MB(d,n1,&B[i*n1*d],false);
  1010. MB=M*MA;
  1011. }
  1012. }
  1013. template <class T>
  1014. void cheb_grad(T* A, int deg, T* B){
  1015. int dim=3;
  1016. int d=deg+1;
  1017. int n1 =(d*(d+1)*(d+2))/6;
  1018. int n1_=(int)(pow((T)d,dim)+0.5);
  1019. Vector<T> A_(n1_); A_.SetZero();
  1020. Vector<T> B_(n1_); B_.SetZero();
  1021. {// Rearrange data
  1022. int indx=0;
  1023. for(int i=0;i<d;i++)
  1024. for(int j=0;i+j<d;j++)
  1025. for(int k=0;i+j+k<d;k++){
  1026. A_[k+(j+i*d)*d]=A[indx];
  1027. indx++;
  1028. }
  1029. }
  1030. for(int l=0;l<dim;l++){
  1031. cheb_diff(&A_[0],d-1,dim,l,&B_[0]);
  1032. {// Rearrange data
  1033. int indx=l*n1;
  1034. for(int i=0;i<d;i++)
  1035. for(int j=0;i+j<d;j++)
  1036. for(int k=0;i+j+k<d;k++){
  1037. B[indx]=B_[k+(j+i*d)*d];
  1038. indx++;
  1039. }
  1040. }
  1041. }
  1042. }
  1043. template <class T>
  1044. void cheb_div(T* A_, int deg, T* B_){
  1045. int dim=3;
  1046. int d=deg+1;
  1047. int n1 =(int)(pow((T)d,dim)+0.5);
  1048. Vector<T> A(n1*dim); A.SetZero();
  1049. Vector<T> B(n1 ); B.SetZero();
  1050. {// Rearrange data
  1051. int indx=0;
  1052. for(int l=0;l<dim;l++)
  1053. for(int i=0;i<d;i++)
  1054. for(int j=0;i+j<d;j++)
  1055. for(int k=0;i+j+k<d;k++){
  1056. A[k+(j+(i+l*d)*d)*d]=A_[indx];
  1057. indx++;
  1058. }
  1059. }
  1060. Matrix<T> MB(n1,1,&B[0],false);
  1061. Matrix<T> MC(n1,1);
  1062. for(int i=0;i<3;i++){
  1063. cheb_diff(&A[n1*i],d-1,3,i,MC[0]);
  1064. MB+=MC;
  1065. }
  1066. {// Rearrange data
  1067. int indx=0;
  1068. for(int i=0;i<d;i++)
  1069. for(int j=0;i+j<d;j++)
  1070. for(int k=0;i+j+k<d;k++){
  1071. B_[indx]=B[k+(j+i*d)*d];
  1072. indx++;
  1073. }
  1074. }
  1075. }
  1076. template <class T>
  1077. void cheb_curl(T* A_, int deg, T* B_){
  1078. int dim=3;
  1079. int d=deg+1;
  1080. int n1 =(int)(pow((T)d,dim)+0.5);
  1081. Vector<T> A(n1*dim); A.SetZero();
  1082. Vector<T> B(n1*dim); B.SetZero();
  1083. {// Rearrange data
  1084. int indx=0;
  1085. for(int l=0;l<dim;l++)
  1086. for(int i=0;i<d;i++)
  1087. for(int j=0;i+j<d;j++)
  1088. for(int k=0;i+j+k<d;k++){
  1089. A[k+(j+(i+l*d)*d)*d]=A_[indx];
  1090. indx++;
  1091. }
  1092. }
  1093. Matrix<T> MC1(n1,1);
  1094. Matrix<T> MC2(n1,1);
  1095. for(int i=0;i<3;i++){
  1096. Matrix<T> MB(n1,1,&B[n1*i],false);
  1097. int j1=(i+1)%3;
  1098. int j2=(i+2)%3;
  1099. cheb_diff(&A[n1*j1],d-1,3,j2,MC1[0]);
  1100. cheb_diff(&A[n1*j2],d-1,3,j1,MC2[0]);
  1101. MB=MC2;
  1102. MB-=MC1;
  1103. }
  1104. {// Rearrange data
  1105. int indx=0;
  1106. for(int l=0;l<dim;l++)
  1107. for(int i=0;i<d;i++)
  1108. for(int j=0;i+j<d;j++)
  1109. for(int k=0;i+j+k<d;k++){
  1110. B_[indx]=B[k+(j+(i+l*d)*d)*d];
  1111. indx++;
  1112. }
  1113. }
  1114. }
  1115. //TODO: Fix number of cheb_coeff to (d+1)*(d+2)*(d+3)/6 for the following functions.
  1116. template <class T>
  1117. void cheb_laplacian(T* A, int deg, T* B){
  1118. int dim=3;
  1119. int d=deg+1;
  1120. int n1=(int)(pow((T)d,dim)+0.5);
  1121. T* C1=new T[n1];
  1122. T* C2=new T[n1];
  1123. Matrix<T> M_(1,n1,C2,false);
  1124. for(int i=0;i<3;i++){
  1125. Matrix<T> M (1,n1,&B[n1*i],false);
  1126. for(int j=0;j<n1;j++) M[0][j]=0;
  1127. for(int j=0;j<3;j++){
  1128. cheb_diff(&A[n1*i],d-1,3,j,C1);
  1129. cheb_diff( C1 ,d-1,3,j,C2);
  1130. M+=M_;
  1131. }
  1132. }
  1133. delete[] C1;
  1134. delete[] C2;
  1135. }
  1136. /*
  1137. * \brief Computes image of the chebyshev interpolation along the specified axis.
  1138. */
  1139. template <class T>
  1140. void cheb_img(T* A, T* B, int deg, int dir, bool neg_){
  1141. int d=deg+1;
  1142. int n1=(int)(pow((T)d,3-dir)+0.5);
  1143. int n2=(int)(pow((T)d, dir)+0.5);
  1144. int indx;
  1145. T sgn,neg;
  1146. neg=(T)(neg_?-1.0:1.0);
  1147. for(int i=0;i<n1;i++){
  1148. indx=i%d;
  1149. sgn=(T)(indx%2?-neg:neg);
  1150. for(int j=0;j<n2;j++){
  1151. B[i*n2+j]=sgn*A[i*n2+j];
  1152. }
  1153. }
  1154. }
  1155. }//end namespace