cheb_utils.txx 32 KB

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