cheb_utils.txx 31 KB

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