/** * \file cheb_node.cpp * \author Dhairya Malhotra, dhairya.malhotra@gmail.com * \date 1-22-2010 * \brief This file contains the implementation of the class Cheb_Node. */ #include #include #include #include namespace pvfmm{ template Cheb_Node::~Cheb_Node(){} template void Cheb_Node::Initialize(TreeNode* parent_, int path2node_, TreeNode::NodeData* data_) { MPI_Node::Initialize(parent_,path2node_,data_); //Set Cheb_Node specific data. NodeData* cheb_data=dynamic_cast(data_); Cheb_Node* parent=dynamic_cast*>(this->Parent()); if(cheb_data!=NULL){ cheb_deg=cheb_data->cheb_deg; input_fn=cheb_data->input_fn; data_dof=cheb_data->data_dof; tol=cheb_data->tol; }else if(parent!=NULL){ cheb_deg=parent->cheb_deg; input_fn=parent->input_fn; data_dof=parent->data_dof; tol=parent->tol; } //Compute Chebyshev approximation. if(this->IsLeaf() && !this->IsGhost()){ if(input_fn!=NULL && data_dof>0){ Real_t s=pow(0.5,this->Depth()); int n1=(int)(pow((Real_t)(cheb_deg+1),this->Dim())+0.5); std::vector coord=cheb_nodes(cheb_deg,this->Dim()); for(int i=0;iCoord()[0]; coord[i*3+1]=coord[i*3+1]*s+this->Coord()[1]; coord[i*3+2]=coord[i*3+2]*s+this->Coord()[2]; } std::vector input_val(n1*data_dof); input_fn(&coord[0],n1,&input_val[0]); Matrix M_val(n1,data_dof,&input_val[0],false); M_val=M_val.Transpose(); cheb_coeff.Resize(((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6*data_dof); cheb_coeff.SetZero(); cheb_approx(&input_val[0], cheb_deg, data_dof, &cheb_coeff[0]); }else if(this->cheb_value.Dim()>0){ size_t n_ptr=this->cheb_coord.Dim()/this->Dim(); assert(n_ptr*data_dof==this->cheb_value.Dim()); points2cheb(cheb_deg,&(this->cheb_coord[0]),&(this->cheb_value[0]), this->cheb_coord.Dim()/this->Dim(),data_dof,this->Coord(), 1.0/(1UL<Depth()), cheb_coeff); } } } template void Cheb_Node::ClearData(){ ChebData().Resize(0); MPI_Node::ClearData(); } template TreeNode* Cheb_Node::NewNode(TreeNode* n_){ Cheb_Node* n=(n_==NULL?new Cheb_Node():static_cast*>(n_)); n->cheb_deg=cheb_deg; n->input_fn=input_fn; n->data_dof=data_dof; n->tol=tol; return MPI_Node::NewNode(n); } template bool Cheb_Node::SubdivCond(){ // Do not subdivide beyond max_depth if(this->Depth()>=this->max_depth-1) return false; if(!this->IsLeaf()){ // If has non-leaf children, then return true. int n=(1UL<Dim()); for(int i=0;i* ch=static_cast*>(this->Child(i)); assert(ch!=NULL); //This should never happen if(!ch->IsLeaf() || ch->IsGhost()) return true; } } else{ // Do not refine ghost leaf nodes. if(this->IsGhost()) return false; } if(MPI_Node::SubdivCond()) return true; if(!this->IsLeaf()){ std::vector val((cheb_deg+1)*(cheb_deg+1)*(cheb_deg+1)*data_dof); std::vector x=cheb_nodes(cheb_deg,1); std::vector y=x; std::vector z=x; Real_t s=pow(0.5,this->Depth()); Real_t* coord=this->Coord(); for(size_t i=0;i tmp_coeff(data_dof*((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6); cheb_approx(&val[0],cheb_deg,data_dof,&tmp_coeff[0]); return (cheb_err(&(tmp_coeff[0]),cheb_deg,data_dof)>tol); }else{ assert(cheb_coeff.Dim()>0); Real_t err=cheb_err(&(cheb_coeff[0]),cheb_deg,data_dof); return (err>tol); } } template void Cheb_Node::Subdivide() { if(!this->IsLeaf()) return; MPI_Node::Subdivide(); if(cheb_deg<0 || cheb_coeff.Dim()==0 || input_fn!=NULL) return; std::vector x(cheb_deg+1); std::vector y(cheb_deg+1); std::vector z(cheb_deg+1); Vector cheb_node=cheb_nodes(cheb_deg,1); Vector val((size_t)pow(cheb_deg+1,this->Dim())*data_dof); Vector child_cheb_coeff[8]; int n=(1UL<Dim()); for(int i=0;iDim())*data_dof); child_cheb_coeff[i].Resize(data_dof*((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6); cheb_approx(&val[0],cheb_deg,data_dof,&(child_cheb_coeff[i][0])); Cheb_Node* child=static_cast*>(this->Child(i)); child->cheb_coeff=child_cheb_coeff[i]; assert(child->cheb_deg==cheb_deg); assert(child->input_fn==input_fn); assert(child->tol==tol); } } template void Cheb_Node::Truncate() { if(cheb_deg<0 || this->IsLeaf()) return MPI_Node::Truncate(); std::vector cheb_node=cheb_nodes(cheb_deg,1); std::vector x=cheb_node; std::vector y=cheb_node; std::vector z=cheb_node; std::vector val((cheb_deg+1)*(cheb_deg+1)*(cheb_deg+1)*data_dof); Real_t s=pow(0.5,this->Depth()); Real_t* coord=this->Coord(); for(size_t i=0;i(&val[0],cheb_deg,data_dof,&cheb_coeff[0]); MPI_Node::Truncate(); } template PackedData Cheb_Node::Pack(bool ghost, void* buff_ptr, size_t offset) { return MPI_Node::Pack(ghost,buff_ptr, offset); } template void Cheb_Node::Unpack(PackedData p0, bool own_data) { MPI_Node::Unpack(p0, own_data); } template void Cheb_Node::VTU_Data(std::vector& coord, std::vector& value, std::vector& connect, std::vector& offset, std::vector& types, int lod){ //for(size_t i=0;ipt_coord.Dim();i++) coord.push_back(this->pt_coord[i]); //for(size_t i=0;ipt_value.Dim();i++) value.push_back(this->pt_value[i]); Real_t* c=this->Coord(); Real_t s=pow(0.5,this->Depth()); size_t point_cnt=coord.size()/COORD_DIM; int gridpt_cnt; std::vector grid_pts; { int cheb_deg=lod; std::vector cheb_node_=cheb_nodes(cheb_deg,1); gridpt_cnt=cheb_node_.size()+2; grid_pts.resize(gridpt_cnt); for(size_t i=0;i x(gridpt_cnt); std::vector y(gridpt_cnt); std::vector z(gridpt_cnt); for(int i=0;iReadVal(x, y, z, val_ptr); //Rearrrange data //(x1,x2,x3,...,y1,y2,...z1,...) => (x1,y1,z1,x2,y2,z2,...) Matrix M(data_dof,gridpt_cnt*gridpt_cnt*gridpt_cnt,val_ptr,false); M=M.Transpose(); } } template void Cheb_Node::Gradient(){ int dim=3;//this->Dim(); if(this->IsLeaf() && ChebData().Dim()>0){ int n3=((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6; Vector coeff(ChebData().Dim()*dim); for(int i=0;idepth); for(size_t i=0;i void Cheb_Node::Divergence(){ int dim=3;//this->Dim(); if(this->IsLeaf() && ChebData().Dim()>0){ assert(data_dof%3==0); int n3=((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6; Vector coeff(ChebData().Dim()/dim); for(int i=0;idepth); for(size_t i=0;i void Cheb_Node::Curl(){ int dim=3;//this->Dim(); if(this->IsLeaf() && ChebData().Dim()>0){ assert(data_dof%dim==0); int n3=((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6; Vector coeff(ChebData().Dim()); for(int i=0;idepth); for(size_t i=0;i void Cheb_Node::read_val(std::vector x,std::vector y, std::vector z, int nx, int ny, int nz, Real_t* val, bool show_ghost/*=true*/){ if(cheb_deg<0) return; Real_t s=0.5*pow(0.5,this->Depth()); Real_t s_inv=1/s; if(this->IsLeaf()){ if(cheb_coeff.Dim()!=(size_t)((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6*data_dof || (this->IsGhost() && !show_ghost)) return; Vector out; std::vector x_=x; std::vector y_=y; std::vector z_=z; for(size_t i=0;iCoord()[0])*s_inv-1.0; for(size_t i=0;iCoord()[1])*s_inv-1.0; for(size_t i=0;iCoord()[2])*s_inv-1.0; cheb_eval(cheb_coeff, cheb_deg, x_, y_, z_, out); for(int l=0;lCoord()[0]+s, this->Coord()[1]+s, this->Coord()[2]+s}; Real_t* indx[3]={&(std::lower_bound(x.begin(),x.end(),coord_[0])[0]), &(std::lower_bound(y.begin(),y.end(),coord_[1])[0]), &(std::lower_bound(z.begin(),z.end(),coord_[2])[0])}; std::vector x1[2]={std::vector(&(x.begin()[0]),indx[0]),std::vector(indx[0],&(x.end()[0]))}; std::vector y1[2]={std::vector(&(y.begin()[0]),indx[1]),std::vector(indx[1],&(y.end()[0]))}; std::vector z1[2]={std::vector(&(z.begin()[0]),indx[2]),std::vector(indx[2],&(z.end()[0]))}; for(int i=0;i<8;i++){ std::vector& x1_=x1[i%2]; std::vector& y1_=y1[(i/2)%2]; std::vector& z1_=z1[(i/4)%2]; if(x1_.size()>0 && y1_.size()>0 && z1_.size()>0){ static_cast*>(this->Child(i))->read_val(x1_,y1_,z1_,nx,ny,nz,&val[(i%2?indx[0]-&x[0]:0)+((i/2)%2?indx[1]-&y[0]:0)*nx+((i/4)%2?indx[2]-&z[0]:0)*nx*ny],show_ghost); } } } }//end namespace