Skip to content
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
14 changes: 8 additions & 6 deletions src/invert/laplace/impls/cyclic/cyclic_laplace.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -38,8 +38,10 @@
#if not BOUT_USE_METRIC_3D

#include "cyclic_laplace.hxx"
#include "bout/assert.hxx"
#include "bout/bout_types.hxx"

#include <bout/assert.hxx>
#include <bout/bout_types.hxx>
#include <bout/boutcomm.hxx>
#include <bout/boutexception.hxx>
#include <bout/constants.hxx>
#include <bout/fft.hxx>
Expand Down Expand Up @@ -591,21 +593,21 @@ void LaplaceCyclic ::verify_solution(const Matrix<dcomplex>& a_ver,
}

if (xproc > 0) {
MPI_Irecv(&rbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 901, MPI_COMM_WORLD,
MPI_Irecv(&rbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 901, BoutComm::get(),
&request[1]);
for (int kz = 0; kz < nsys; kz++) {
sbufdown[kz] = x_ver(kz, 1);
}
MPI_Isend(&sbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 900, MPI_COMM_WORLD,
MPI_Isend(&sbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 900, BoutComm::get(),
&request[0]);
}
if (xproc < nprocs - 1) {
MPI_Irecv(&rbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 900, MPI_COMM_WORLD,
MPI_Irecv(&rbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 900, BoutComm::get(),
&request[3]);
for (int kz = 0; kz < nsys; kz++) {
sbufup[kz] = x_ver(kz, nx);
}
MPI_Isend(&sbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 901, MPI_COMM_WORLD,
MPI_Isend(&sbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 901, BoutComm::get(),
&request[2]);
}

Expand Down
8 changes: 4 additions & 4 deletions src/invert/laplace/impls/pcr/pcr.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -1041,21 +1041,21 @@ void LaplacePCR ::verify_solution(const Matrix<dcomplex>& a_ver,
}

if (xproc > 0) {
MPI_Irecv(&rbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 901, MPI_COMM_WORLD,
MPI_Irecv(&rbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 901, BoutComm::get(),
&request[1]);
for (int kz = 0; kz < nsys; kz++) {
sbufdown[kz] = x_ver(kz, 1);
}
MPI_Isend(&sbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 900, MPI_COMM_WORLD,
MPI_Isend(&sbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 900, BoutComm::get(),
&request[0]);
}
if (xproc < nprocs - 1) {
MPI_Irecv(&rbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 900, MPI_COMM_WORLD,
MPI_Irecv(&rbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 900, BoutComm::get(),
&request[3]);
for (int kz = 0; kz < nsys; kz++) {
sbufup[kz] = x_ver(kz, nx);
}
MPI_Isend(&sbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 901, MPI_COMM_WORLD,
MPI_Isend(&sbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 901, BoutComm::get(),
&request[2]);
}

Expand Down
8 changes: 4 additions & 4 deletions src/invert/laplace/impls/pcr_thomas/pcr_thomas.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -1140,21 +1140,21 @@ void LaplacePCR_THOMAS ::verify_solution(const Matrix<dcomplex>& a_ver,
}

if (xproc > 0) {
MPI_Irecv(&rbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 901, MPI_COMM_WORLD,
MPI_Irecv(&rbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 901, BoutComm::get(),
&request[1]);
for (int kz = 0; kz < nsys; kz++) {
sbufdown[kz] = x_ver(kz, 1);
}
MPI_Isend(&sbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 900, MPI_COMM_WORLD,
MPI_Isend(&sbufdown[0], nsys, MPI_DOUBLE_COMPLEX, myrank - 1, 900, BoutComm::get(),
&request[0]);
}
if (xproc < nprocs - 1) {
MPI_Irecv(&rbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 900, MPI_COMM_WORLD,
MPI_Irecv(&rbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 900, BoutComm::get(),
&request[3]);
for (int kz = 0; kz < nsys; kz++) {
sbufup[kz] = x_ver(kz, nx);
}
MPI_Isend(&sbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 901, MPI_COMM_WORLD,
MPI_Isend(&sbufup[0], nsys, MPI_DOUBLE_COMPLEX, myrank + 1, 901, BoutComm::get(),
&request[2]);
}

Expand Down
10 changes: 2 additions & 8 deletions src/mesh/interpolation/hermite_spline_xz.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,7 @@
**************************************************************************/

#include "../impls/bout/boutmesh.hxx"
#include "bout/bout.hxx"
#include "bout/globals.hxx"
#include "bout/index_derivs_interface.hxx"
#include "bout/interpolation_xz.hxx"
Expand Down Expand Up @@ -133,16 +134,9 @@ XZHermiteSpline::XZHermiteSpline(int y_offset, Mesh* mesh)
#ifdef HS_USE_PETSC
petsclib = new PetscLib(
&Options::root()["mesh:paralleltransform:xzinterpolation:hermitespline"]);
// MatCreate(MPI_Comm comm,Mat *A)
// MatCreate(MPI_COMM_WORLD, &petscWeights);
// MatSetSizes(petscWeights, m, m, M, M);
// PetscErrorCode MatCreateAIJ(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M,
// PetscInt N, PetscInt d_nz, const PetscInt d_nnz[],
// PetscInt o_nz, const PetscInt o_nnz[], Mat *A)
// MatSetSizes(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)
const int m = localmesh->LocalNx * localmesh->LocalNy * localmesh->LocalNz;
const int M = m * localmesh->getNXPE() * localmesh->getNYPE();
MatCreateAIJ(MPI_COMM_WORLD, m, m, M, M, 16, nullptr, 16, nullptr, &petscWeights);
MatCreateAIJ(BoutComm::get(), m, m, M, M, 16, nullptr, 16, nullptr, &petscWeights);
Copy link
Copy Markdown
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

warning: no header providing "MatCreateAIJ" is directly included [misc-include-cleaner]

src/mesh/interpolation/hermite_spline_xz.cxx:28:

- #include <vector>
+ #include <petscmat.h>
+ #include <vector>

#endif
#endif
#ifndef HS_USE_PETSC
Expand Down
Loading