Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Multipole: only integrate on root rank #29

Merged
merged 1 commit into from
Jun 28, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 2 additions & 1 deletion Multipole/src/multipole.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -155,7 +155,8 @@ extern "C" void Multipole_Calc(CCTK_ARGUMENTS) {
// Intergate of conj(sYlm)*F*sin(theta) over the sphere at radiusr[i]
for (int l = 0; l <= l_max; l++) {
for (int m = -l; m <= l; m++) {
g_sphere->integrate(g_sphere->getRealY()[si][l][m + l],
g_sphere->integrate(CCTK_PASS_CTOC,
g_sphere->getRealY()[si][l][m + l],
g_sphere->getImagY()[si][l][m + l],
g_sphere->getRealF(), g_sphere->getImagF(),
&modes(v, i, l, m, 0), &modes(v, i, l, m, 1));
Expand Down
73 changes: 39 additions & 34 deletions Multipole/src/surface.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -66,50 +66,55 @@ void Surface::interpolate(CCTK_ARGUMENTS, int realFieldIndex,
}

// Take the integral of conj(array1)*array2*sin(th)
void Surface::integrate(const std::vector<CCTK_REAL> &array1r,
void Surface::integrate(CCTK_ARGUMENTS, const std::vector<CCTK_REAL> &array1r,
const std::vector<CCTK_REAL> &array1i,
const std::vector<CCTK_REAL> &array2r,
const std::vector<CCTK_REAL> &array2i, CCTK_REAL *outRe,
CCTK_REAL *outIm) {
DECLARE_CCTK_PARAMETERS;

std::vector<CCTK_REAL> fReal(theta_.size());
std::vector<CCTK_REAL> fImag(theta_.size());
if (CCTK_MyProc(cctkGH) == 0) {

// integrand: conj(array1)*array2*sin(th)
for (size_t i = 0; i < theta_.size(); ++i) {
fReal[i] = (array1r[i] * array2r[i] + array1i[i] * array2i[i]) *
std::sin(theta_[i]);
fImag[i] = (array1r[i] * array2i[i] - array1i[i] * array2r[i]) *
std::sin(theta_[i]);
}
std::vector<CCTK_REAL> fReal(theta_.size());
std::vector<CCTK_REAL> fImag(theta_.size());

if (CCTK_Equals(integration_method, "midpoint")) {
*outRe = Midpoint2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm = Midpoint2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else if (CCTK_Equals(integration_method, "trapezoidal")) {
*outRe =
Trapezoidal2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm =
Trapezoidal2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else if (CCTK_Equals(integration_method, "Simpson")) {
if (nPhi_ % 2 != 0 || nTheta_ % 2 != 0) {
CCTK_WARN(CCTK_WARN_ABORT, "The Simpson integration method requires even "
"nTheta_ and even nPhi_");
// integrand: conj(array1)*array2*sin(th)
for (size_t i = 0; i < theta_.size(); ++i) {
fReal[i] = (array1r[i] * array2r[i] + array1i[i] * array2i[i]) *
std::sin(theta_[i]);
fImag[i] = (array1r[i] * array2i[i] - array1i[i] * array2r[i]) *
std::sin(theta_[i]);
}
*outRe = Simpson2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm = Simpson2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else if (CCTK_Equals(integration_method, "DriscollHealy")) {
if (nTheta_ % 2 != 0) {
CCTK_WARN(CCTK_WARN_ABORT,
"The Driscoll&Healy integration method requires even nTheta_");

if (CCTK_Equals(integration_method, "midpoint")) {
*outRe = Midpoint2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm = Midpoint2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else if (CCTK_Equals(integration_method, "trapezoidal")) {
*outRe =
Trapezoidal2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm =
Trapezoidal2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else if (CCTK_Equals(integration_method, "Simpson")) {
if (nPhi_ % 2 != 0 || nTheta_ % 2 != 0) {
CCTK_WARN(CCTK_WARN_ABORT,
"The Simpson integration method requires even "
"nTheta_ and even nPhi_");
}
*outRe = Simpson2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm = Simpson2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else if (CCTK_Equals(integration_method, "DriscollHealy")) {
if (nTheta_ % 2 != 0) {
CCTK_WARN(
CCTK_WARN_ABORT,
"The Driscoll&Healy integration method requires even nTheta_");
}
*outRe =
DriscollHealy2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm =
DriscollHealy2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else {
CCTK_WARN(CCTK_WARN_ABORT, "internal error");
}
*outRe =
DriscollHealy2DIntegral(fReal.data(), nTheta_, nPhi_, dTheta_, dPhi_);
*outIm =
DriscollHealy2DIntegral(fImag.data(), nTheta_, nPhi_, dTheta_, dPhi_);
} else {
CCTK_WARN(CCTK_WARN_ABORT, "internal error");
}
}

Expand Down
2 changes: 1 addition & 1 deletion Multipole/src/surface.hxx
Original file line number Diff line number Diff line change
Expand Up @@ -60,7 +60,7 @@ public:
void interpolate(CCTK_ARGUMENTS, int realFieldIndex, int imagFieldIndex);

// Take the integral of conj(array1)*array2*sin(th)
void integrate(const std::vector<CCTK_REAL> &array1r,
void integrate(CCTK_ARGUMENTS, const std::vector<CCTK_REAL> &array1r,
const std::vector<CCTK_REAL> &array1i,
const std::vector<CCTK_REAL> &array2r,
const std::vector<CCTK_REAL> &array2i, CCTK_REAL *outre,
Expand Down
Loading