Skip to content

Commit

Permalink
Merge remote-tracking branch 'origin/develop' into cuda-pts
Browse files Browse the repository at this point in the history
  • Loading branch information
blackwer committed Dec 6, 2022
2 parents a2238ec + 5ac4e71 commit 1c7af14
Show file tree
Hide file tree
Showing 4 changed files with 26 additions and 11 deletions.
8 changes: 6 additions & 2 deletions examples/src/example2-f.f90
Original file line number Diff line number Diff line change
Expand Up @@ -31,7 +31,9 @@ program main
subroutine fn_input(coord, n, val)
implicit none
integer*8 :: n, i, dof
real*8 :: coord(n*3), val(n*3), L, r_2
real*8, intent(in) :: coord(n*3)
real*8, intent(out) :: val(n*3)
real*8 :: L, r_2

dof = 3
L = 125
Expand All @@ -46,7 +48,9 @@ subroutine fn_input(coord, n, val)
subroutine fn_poten(coord, n, val)
implicit none
integer*8 :: n, i, dof
real*8 :: coord(n*3), val(n*3), L, r_2
real*8, intent(in) :: coord(n*3)
real*8, intent(out) :: val(n*3)
real*8 :: L, r_2

dof = 3
L = 125
Expand Down
6 changes: 3 additions & 3 deletions include/fmm_pts.txx
Original file line number Diff line number Diff line change
Expand Up @@ -1073,9 +1073,6 @@ Matrix<typename FMMNode::Real_t>& FMM_Pts<FMMNode>::Precomp(int level, Mat_Type

int xlow,xhigh,ylow,yhigh,zlow,zhigh;
switch(mat_indx){
case BoundaryType::FreeSpace :
xlow=0;xhigh=0;ylow=0;yhigh=0;zlow=0;zhigh=0;
break;
case BoundaryType::PX :
xlow=-2;xhigh=2;
ylow=0;yhigh=0;
Expand All @@ -1091,6 +1088,9 @@ Matrix<typename FMMNode::Real_t>& FMM_Pts<FMMNode>::Precomp(int level, Mat_Type
ylow=-2;yhigh=2;
zlow=-2;zhigh=2;
break;
default:
xlow=0;xhigh=0;ylow=0;yhigh=0;zlow=0;zhigh=0;
break;
}

for(int x0=xlow;x0<=xhigh;x0++)
Expand Down
5 changes: 5 additions & 0 deletions include/mpi_tree.txx
Original file line number Diff line number Diff line change
Expand Up @@ -1165,6 +1165,11 @@ void MPI_Tree<TreeNode>::SetColleagues(BoundaryType bndry, Node_t* node){
ylow=-1;yhigh=+1;
zlow=0;zhigh=0;
break;
default:
xlow=0;xhigh=0;
ylow=0;yhigh=0;
zlow=0;zhigh=0;
break;
}

for(long i0=xlow;i0<=xhigh;i0++)
Expand Down
18 changes: 12 additions & 6 deletions src/pvfmm-wrapper.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,12 @@
extern "C" { // Volume FM
#endif

#ifdef PVFMM_EXTENDED_BC
#define PVFMM_FULL_PERIODIC pvfmm::PXYZ
#else
#define PVFMM_FULL_PERIODIC pvfmm::Periodic
#endif

void* PVFMMCreateVolumeFMMF(int m, int q, enum PVFMMKernel kernel, MPI_Comm comm) {
const pvfmm::Kernel<float>* ker = nullptr;
if (kernel == PVFMMLaplacePotential ) ker = &pvfmm::LaplaceKernel <float>::potential();
Expand All @@ -30,7 +36,7 @@ void* PVFMMCreateVolumeTreeF(int cheb_deg, int data_dim, void (*fn_ptr)(const fl
fn_ptr(coord, n, out, fn_ctx);
};

auto* tree = ChebFMM_CreateTree(cheb_deg, data_dim, fn_ptr_, trg_coord_, comm, tol, max_pts, periodic?pvfmm::Periodic:pvfmm::FreeSpace, init_depth);
auto* tree = ChebFMM_CreateTree(cheb_deg, data_dim, fn_ptr_, trg_coord_, comm, tol, max_pts, periodic?PVFMM_FULL_PERIODIC:pvfmm::FreeSpace, init_depth);
//tree->Write2File("vis",4);

return (void*)tree;
Expand All @@ -46,7 +52,7 @@ void* PVFMMCreateVolumeTreeFromCoeffF(long n_nodes, int cheb_deg, int data_dim,
#pragma omp parallel for schedule(static)
for (long i = 0; i < (long)trg_coord_.size(); i++) trg_coord_[i] = trg_coord[i];

auto* tree = ChebFMM_CreateTree(cheb_deg, node_coord_, fn_coeff_, trg_coord_, comm, periodic?pvfmm::Periodic:pvfmm::FreeSpace);
auto* tree = ChebFMM_CreateTree(cheb_deg, node_coord_, fn_coeff_, trg_coord_, comm, periodic?PVFMM_FULL_PERIODIC:pvfmm::FreeSpace);
//tree->Write2File("vis",4);

return (void*)tree;
Expand Down Expand Up @@ -138,7 +144,7 @@ void* PVFMMCreateVolumeTreeD(int cheb_deg, int data_dim, void (*fn_ptr)(const do
fn_ptr(coord, n, out, fn_ctx);
};

auto* tree = ChebFMM_CreateTree(cheb_deg, data_dim, fn_ptr_, trg_coord_, comm, tol, max_pts, periodic?pvfmm::Periodic:pvfmm::FreeSpace, init_depth);
auto* tree = ChebFMM_CreateTree(cheb_deg, data_dim, fn_ptr_, trg_coord_, comm, tol, max_pts, periodic?PVFMM_FULL_PERIODIC:pvfmm::FreeSpace, init_depth);
//tree->Write2File("vis",4);

return (void*)tree;
Expand All @@ -154,7 +160,7 @@ void* PVFMMCreateVolumeTreeFromCoeffD(long n_nodes, int cheb_deg, int data_dim,
#pragma omp parallel for schedule(static)
for (long i = 0; i < (long)trg_coord_.size(); i++) trg_coord_[i] = trg_coord[i];

auto* tree = ChebFMM_CreateTree(cheb_deg, node_coord_, fn_coeff_, trg_coord_, comm, periodic?pvfmm::Periodic:pvfmm::FreeSpace);
auto* tree = ChebFMM_CreateTree(cheb_deg, node_coord_, fn_coeff_, trg_coord_, comm, periodic?PVFMM_FULL_PERIODIC:pvfmm::FreeSpace);
//tree->Write2File("vis",4);

return (void*)tree;
Expand Down Expand Up @@ -242,7 +248,7 @@ void pvfmmcreatevolumetreef_(void** ctx, const int32_t* cheb_deg, const int32_t*
};

const MPI_Comm comm = MPI_Comm_f2c(*fcomm);
auto* tree = ChebFMM_CreateTree(*cheb_deg, *data_dim, fn_ptr_, trg_coord_, comm, *tol, *max_pts, (*periodic)==0?pvfmm::FreeSpace:pvfmm::Periodic, *init_depth);
auto* tree = ChebFMM_CreateTree(*cheb_deg, *data_dim, fn_ptr_, trg_coord_, comm, *tol, *max_pts, (*periodic)==0?pvfmm::FreeSpace:PVFMM_FULL_PERIODIC, *init_depth);
//tree->Write2File("vis",4);

(*ctx) = (void*)tree;
Expand Down Expand Up @@ -304,7 +310,7 @@ void pvfmmcreatevolumetreed_(void** ctx, const int32_t* cheb_deg, const int32_t*
};

const MPI_Comm comm = MPI_Comm_f2c(*fcomm);
auto* tree = ChebFMM_CreateTree(*cheb_deg, *data_dim, fn_ptr_, trg_coord_, comm, *tol, *max_pts, (*periodic)==0?pvfmm::FreeSpace:pvfmm::Periodic, *init_depth);
auto* tree = ChebFMM_CreateTree(*cheb_deg, *data_dim, fn_ptr_, trg_coord_, comm, *tol, *max_pts, (*periodic)==0?pvfmm::FreeSpace:PVFMM_FULL_PERIODIC, *init_depth);
//tree->Write2File("vis",4);

(*ctx) = (void*)tree;
Expand Down

0 comments on commit 1c7af14

Please sign in to comment.