cheb_utils.txx 27 KB

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