cheb_utils.txx 35 KB

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