Skip to content

Commit

Permalink
AztecOO: Implement setter/getter for new global_comm wrapper
Browse files Browse the repository at this point in the history
  • Loading branch information
JacobDomagala committed Sep 1, 2023
1 parent a37e42d commit db28a3b
Show file tree
Hide file tree
Showing 6 changed files with 155 additions and 25 deletions.
62 changes: 62 additions & 0 deletions packages/aztecoo/src/AztecOO_GlobalComm.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,62 @@

/*
//@HEADER
// ***********************************************************************
//
// AztecOO: An Object-Oriented Aztec Linear Solver Package
// Copyright (2002) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ***********************************************************************
//@HEADER
*/

#include "AztecOO_GlobalComm.h"

pthread_mutex_t mpi_mutex = PTHREAD_MUTEX_INITIALIZER;
MPI_Comm Global_AztecOO_Comm = MPI_COMM_WORLD;

void AZ_initialize_global_comm(MPI_Comm comm) {
pthread_mutex_lock(&mpi_mutex);
Global_AztecOO_Comm = comm;
pthread_mutex_unlock(&mpi_mutex);
}

MPI_Comm AZ_get_global_comm() {
MPI_Comm tmp;
pthread_mutex_lock(&mpi_mutex);
tmp = Global_AztecOO_Comm;
pthread_mutex_unlock(&mpi_mutex);
return tmp;
}
69 changes: 69 additions & 0 deletions packages/aztecoo/src/AztecOO_GlobalComm.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,69 @@

/*
//@HEADER
// ***********************************************************************
//
// AztecOO: An Object-Oriented Aztec Linear Solver Package
// Copyright (2002) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ***********************************************************************
//@HEADER
*/

#ifndef AZTECOO_GLOBAL_COMM_HPP
#define AZTECOO_GLOBAL_COMM_HPP

#include "AztecOO_ConfigDefs.h"

#ifdef HAVE_MPI

#include <mpi.h>
#include <pthread.h>

#ifdef __cplusplus
extern "C" {
#endif

// Function prototypes
void AZ_initialize_global_comm(MPI_Comm comm);
MPI_Comm AZ_get_global_comm();

#ifdef __cplusplus
} // extern "C"
#endif

#endif // HAVE_MPI

#endif // AZTECOO_GLOBAL_COMM_HPP
12 changes: 7 additions & 5 deletions packages/aztecoo/src/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -50,7 +50,7 @@ SET(HEADERS ${HEADERS}
)

SET(SOURCES ${SOURCES}
az_comm.c
az_comm.c
)

#
Expand All @@ -73,6 +73,7 @@ SET(HEADERS ${HEADERS}
AztecOO_StatusTestResNorm.h
AztecOO_ConditionNumber.h
AztecOO_ConfigDefs.h
AztecOO_GlobalComm.h
)

SET(SOURCES ${SOURCES}
Expand All @@ -88,10 +89,11 @@ SET(SOURCES ${SOURCES}
AztecOO_StatusTestMaxIters.cpp
AztecOO_StatusTestResNorm.cpp
AztecOO_ConditionNumber.cpp
AztecOO_GlobalComm.cpp
)

SET(HEADERS ${HEADERS}
AztecOO_Version.h
AztecOO_Version.h
)

SET(SOURCES ${SOURCES}
Expand All @@ -109,12 +111,12 @@ IF (TPL_ENABLE_MPI)

SET(SOURCES ${SOURCES}
md_timer_mpi.c
md_wrap_mpi_c.c
md_wrap_mpi_c.c
)
ELSE()
SET(HEADERS ${HEADERS}
md_timer_generic.c
md_wrap_scalar_c.c
md_wrap_scalar_c.c
)
ENDIF()

Expand Down Expand Up @@ -151,7 +153,7 @@ SET(SOURCES ${SOURCES}
az_precond.c
az_scaling.c
az_solve.c
az_sort.c
az_sort.c
)

IF (${PROJECT_NAME}_ENABLE_Fortran)
Expand Down
2 changes: 0 additions & 2 deletions packages/aztecoo/src/az_comm.c
Original file line number Diff line number Diff line change
Expand Up @@ -81,8 +81,6 @@

extern int az_iterate_id;

// Replaces the occurences of hardcoded MPI_COMM_WORLD
MPI_AZComm az_global_mpi_comm = AZ_MPI_COMM_WORLD;

/******************************************************************************/
/******************************************************************************/
Expand Down
7 changes: 3 additions & 4 deletions packages/aztecoo/src/az_old_matvec_mult.c
Original file line number Diff line number Diff line change
Expand Up @@ -44,8 +44,7 @@
#include <stdio.h>
#include <float.h>
#include "az_aztec.h"

extern MPI_AZComm az_global_mpi_comm;
#include "AztecOO_GlobalComm.h"

void AZ_matvec_mult(double *val, int *indx, int *bindx, int *rpntr, int *cpntr,
int *bpntr, double *b, register double *c,
Expand Down Expand Up @@ -122,9 +121,9 @@ void AZ_matvec_mult(double *val, int *indx, int *bindx, int *rpntr, int *cpntr,
Amat.aux_matrix = NULL;
Amat.matrix_type = data_org[AZ_matrix_type];
#ifdef AZTEC_MPI
AZ_set_comm(proc_config, az_global_mpi_comm);
AZ_set_comm(proc_config, AZ_get_global_comm());
if (first_time == 1) {
AZ_set_proc_config(proc_config, az_global_mpi_comm);
AZ_set_proc_config(proc_config, AZ_get_global_comm());
#else
if (first_time == 1) {
AZ_set_proc_config(proc_config, AZ_NOT_MPI);
Expand Down
28 changes: 14 additions & 14 deletions packages/aztecoo/src/md_wrap_mpi_c.c
Original file line number Diff line number Diff line change
Expand Up @@ -60,24 +60,24 @@
#include <stdio.h>
#include <mpi.h>

#include "AztecOO_GlobalComm.h"

int gl_rbuf = 3;
int gl_sbuf = 3;
/******************************************************************************/
/******************************************************************************/
/******************************************************************************/
int the_proc_name = -1;

// declared in az_comm.c (we use MPI_Comm instead of MPI_AZComm as it's guaranteed AztecOO is built with MPI)
extern MPI_Comm az_global_mpi_comm;

void get_parallel_info(int *proc, int *nprocs, int *dim)

{

/* local variables */

MPI_Comm_size(az_global_mpi_comm, nprocs);
MPI_Comm_rank(az_global_mpi_comm, proc);
MPI_Comm_size(AZ_get_global_comm(), nprocs);
MPI_Comm_rank(AZ_get_global_comm(), proc);
*dim = 0;
the_proc_name = *proc;

Expand All @@ -102,11 +102,11 @@ int md_read(char *buf, int bytes, int *source, int *type, int *flag)
if (*source == -1) *source = MPI_ANY_SOURCE;

if (bytes == 0) {
err = MPI_Recv(&gl_rbuf, 1, MPI_BYTE, *source, *type, az_global_mpi_comm,
err = MPI_Recv(&gl_rbuf, 1, MPI_BYTE, *source, *type, AZ_get_global_comm(),
&status);
}
else {
err = MPI_Recv(buf, bytes, MPI_BYTE, *source, *type, az_global_mpi_comm,
err = MPI_Recv(buf, bytes, MPI_BYTE, *source, *type, AZ_get_global_comm(),
&status);
}

Expand All @@ -132,10 +132,10 @@ int md_write(char *buf, int bytes, int dest, int type, int *flag)
int err;

if (bytes == 0) {
err = MPI_Send(&gl_sbuf, 1, MPI_BYTE, dest, type, az_global_mpi_comm);
err = MPI_Send(&gl_sbuf, 1, MPI_BYTE, dest, type, AZ_get_global_comm());
}
else {
err = MPI_Send(buf, bytes, MPI_BYTE, dest, type, az_global_mpi_comm);
err = MPI_Send(buf, bytes, MPI_BYTE, dest, type, AZ_get_global_comm());
}

if (err != 0) (void) fprintf(stderr, "MPI_Send error = %d\n", err);
Expand Down Expand Up @@ -185,11 +185,11 @@ int md_wrap_iread(void *buf, int bytes, int *source, int *type,
if (*source == -1) *source = MPI_ANY_SOURCE;

if (bytes == 0) {
err = MPI_Irecv(&gl_rbuf, 1, MPI_BYTE, *source, *type, az_global_mpi_comm,
err = MPI_Irecv(&gl_rbuf, 1, MPI_BYTE, *source, *type, AZ_get_global_comm(),
request);
}
else {
err = MPI_Irecv(buf, bytes, MPI_BYTE, *source, *type, az_global_mpi_comm,
err = MPI_Irecv(buf, bytes, MPI_BYTE, *source, *type, AZ_get_global_comm(),
request);
}

Expand Down Expand Up @@ -234,10 +234,10 @@ int md_wrap_write(void *buf, int bytes, int dest, int type, int *flag)
int err = 0;

if (bytes == 0) {
err = MPI_Send(&gl_sbuf, 1, MPI_BYTE, dest, type, az_global_mpi_comm);
err = MPI_Send(&gl_sbuf, 1, MPI_BYTE, dest, type, AZ_get_global_comm());
}
else {
err = MPI_Send(buf, bytes, MPI_BYTE, dest, type, az_global_mpi_comm);
err = MPI_Send(buf, bytes, MPI_BYTE, dest, type, AZ_get_global_comm());
}

return err;
Expand Down Expand Up @@ -335,11 +335,11 @@ int md_wrap_iwrite(void *buf, int bytes, int dest, int type, int *flag,
int err = 0;

if (bytes == 0) {
err = MPI_Isend(&gl_sbuf, 1, MPI_BYTE, dest, type, az_global_mpi_comm,
err = MPI_Isend(&gl_sbuf, 1, MPI_BYTE, dest, type, AZ_get_global_comm(),
request);
}
else {
err = MPI_Isend(buf, bytes, MPI_BYTE, dest, type, az_global_mpi_comm,
err = MPI_Isend(buf, bytes, MPI_BYTE, dest, type, AZ_get_global_comm(),
request);
}

Expand Down

0 comments on commit db28a3b

Please sign in to comment.