Skip to content

Commit d382539

Browse files
committed
Move PETSc Initialize and Finalize to amrex::Initialize and Finalize
Close #5030.
1 parent 55c25f7 commit d382539

File tree

2 files changed

+37
-11
lines changed

2 files changed

+37
-11
lines changed

Src/Base/AMReX.cpp

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -29,6 +29,10 @@
2929
#include <AMReX_Sundials.H>
3030
#endif
3131

32+
#ifdef AMREX_USE_PETSC
33+
#include <petscsys.h>
34+
#endif
35+
3236
#ifdef AMREX_USE_CUPTI
3337
#include <AMReX_CuptiTrace.H>
3438
#endif
@@ -135,6 +139,12 @@ namespace {
135139
std::vector<std::string> command_arguments;
136140
}
137141

142+
#ifdef AMREX_USE_PETSC
143+
namespace {
144+
bool petsc_initialized_by_amrex = false;
145+
}
146+
#endif
147+
138148
namespace {
139149
std::streamsize prev_out_precision;
140150
std::streamsize prev_err_precision;
@@ -762,6 +772,20 @@ amrex::Initialize (int& argc, char**& argv, bool build_parm_parse,
762772
}
763773
#endif
764774

775+
#ifdef AMREX_USE_PETSC
776+
{
777+
PetscBool petsc_is_initialized = PETSC_FALSE;
778+
PetscInitialized(&petsc_is_initialized);
779+
if (!petsc_is_initialized) {
780+
PetscErrorCode ierr = PetscInitialize(&argc, &argv, nullptr, nullptr);
781+
if (ierr != 0) {
782+
amrex::Abort("PetscInitialize failed");
783+
}
784+
petsc_initialized_by_amrex = true;
785+
}
786+
}
787+
#endif
788+
765789
#ifdef AMREX_USE_SUNDIALS
766790
sundials::Initialize(amrex::OpenMP::get_max_threads());
767791
#endif
@@ -797,6 +821,13 @@ amrex::Finalize (amrex::AMReX* pamrex)
797821

798822
AMReX::erase(pamrex);
799823

824+
#ifdef AMREX_USE_PETSC
825+
if (petsc_initialized_by_amrex) {
826+
PetscFinalize();
827+
petsc_initialized_by_amrex = false;
828+
}
829+
#endif
830+
800831
#ifdef AMREX_USE_HYPRE
801832
if (init_hypre) { HYPRE_Finalize(); }
802833
#endif

Src/Extern/PETSc/AMReX_PETSc.cpp

Lines changed: 6 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -87,9 +87,6 @@ PETScABecLap::PETScABecLap (const BoxArray& grids, const DistributionMapping& dm
8787

8888
diaginv.define(grids,dmap,ncomp,0);
8989

90-
PETSC_COMM_WORLD = comm;
91-
PetscInitialize(nullptr, nullptr, nullptr, nullptr);
92-
9390
solver = std::make_unique<amrex_KSP>();
9491
A = std::make_unique<amrex_Mat>();
9592
b = std::make_unique<amrex_Vec>();
@@ -106,8 +103,6 @@ PETScABecLap::~PETScABecLap ()
106103
m_factory = nullptr;
107104
m_bndry = nullptr;
108105
m_maxorder = -1;
109-
110-
PetscFinalize();
111106
}
112107

113108
void
@@ -181,8 +176,8 @@ void
181176
PETScABecLap::prepareSolver ()
182177
{
183178
int num_procs, myid;
184-
MPI_Comm_size(PETSC_COMM_WORLD, &num_procs);
185-
MPI_Comm_rank(PETSC_COMM_WORLD, &myid);
179+
MPI_Comm_size(comm, &num_procs);
180+
MPI_Comm_rank(comm, &myid);
186181

187182
const BoxArray& ba = acoefs.boxArray();
188183
const DistributionMapping& dm = acoefs.DistributionMap();
@@ -341,7 +336,7 @@ PETScABecLap::prepareSolver ()
341336
Vector<PetscInt> ncells_allprocs(num_procs);
342337
MPI_Allgather(&ncells_proc, sizeof(PetscInt), MPI_CHAR,
343338
ncells_allprocs.data(), sizeof(PetscInt), MPI_CHAR,
344-
PETSC_COMM_WORLD);
339+
comm);
345340
PetscInt proc_begin = 0;
346341
for (int i = 0; i < myid; ++i) {
347342
proc_begin += ncells_allprocs[i];
@@ -400,7 +395,7 @@ PETScABecLap::prepareSolver ()
400395
PetscInt d_nz = (eb_stencil_size + regular_stencil_size) / 2;
401396
// estimated amount of block off diag elements
402397
PetscInt o_nz = d_nz / 2;
403-
MatCreate(PETSC_COMM_WORLD, &A->a);
398+
MatCreate(comm, &A->a);
404399
MatSetType(A->a, MATMPIAIJ);
405400
MatSetSizes(A->a, ncells_proc, ncells_proc, ncells_world, ncells_world);
406401
MatMPIAIJSetPreallocation(A->a, d_nz, nullptr, o_nz, nullptr );
@@ -587,7 +582,7 @@ PETScABecLap::prepareSolver ()
587582
MatAssemblyBegin(A->a, MAT_FINAL_ASSEMBLY);
588583
MatAssemblyEnd(A->a, MAT_FINAL_ASSEMBLY);
589584
// create solver
590-
KSPCreate(PETSC_COMM_WORLD, &solver->a);
585+
KSPCreate(comm, &solver->a);
591586
KSPSetOperators(solver->a, A->a, A->a);
592587

593588
// Set up preconditioner
@@ -603,7 +598,7 @@ PETScABecLap::prepareSolver ()
603598

604599
// we are not using command line options KSPSetFromOptions(solver->a);
605600
// create b & x
606-
VecCreateMPI(PETSC_COMM_WORLD, ncells_proc, ncells_world, &x->a);
601+
VecCreateMPI(comm, ncells_proc, ncells_world, &x->a);
607602
VecDuplicate(x->a, &b->a);
608603
}
609604

0 commit comments

Comments
 (0)