Skip to content

Commit

Permalink
build also with clang
Browse files Browse the repository at this point in the history
  • Loading branch information
Thomas-Ulrich committed Sep 11, 2024
1 parent cd642e7 commit 82d2ee5
Show file tree
Hide file tree
Showing 5 changed files with 65 additions and 61 deletions.
87 changes: 45 additions & 42 deletions .github/workflows/build-tandem.yml
Original file line number Diff line number Diff line change
Expand Up @@ -3,62 +3,60 @@ name: tandem-build
on:
push:

jobs:
build:
runs-on: ubuntu-latest
env:
PETSC_CACHE_FOLDER: ${{ github.workspace }}/.petsc
COMPILER: gcc-12
ARCHITECTURE: arch-linux-c-opt
env:
PETSC_CACHE_FOLDER: ${{ github.workspace }}/.petsc
ARCHITECTURE: arch-linux-c-opt

jobs:
build-and-test:
runs-on: ubuntu-24.04
strategy:
fail-fast: false
matrix:
compiler:
- cc: gcc-13
cxx: g++-13
- cc: clang-18
cxx: clang++-18
steps:
- name: Checkout Repository
uses: actions/checkout@v4
with:
submodules: recursive

- name: Get Ubuntu and Petsc Version
- name: Install Dependencies and Compilers
run: |
sudo apt-get update
sudo apt-get install -y \
gcc-13 g++-13 \
cmake \
openmpi-bin libopenmpi-dev \
libmetis-dev libparmetis-dev \
libeigen3-dev python3-numpy \
libopenblas0 libopenblas-dev \
liblua5.3-0 liblua5.3-dev \
libomp-dev libgomp1
# libxsmm
git clone --depth 1 --branch 1.17 https://github.com/libxsmm/libxsmm.git gitbuild && \
mkdir -p gitbuild && cd gitbuild && \
make generator -j $(nproc) && sudo cp bin/libxsmm_gemm_generator /usr/bin && \
cd .. && rm -rf gitbuild
- name: Get Ubuntu and PETSc Version
id: get_version
run: |
UBUNTU_VERSION=$(lsb_release -r | awk '{print $2}')
echo "Ubuntu Version: $UBUNTU_VERSION"
echo "UBUNTU_VERSION=$UBUNTU_VERSION" >> $GITHUB_ENV
petsc_version=$(curl -s https://pypi.org/pypi/petsc/json | jq -r '.info.version')
echo "PETSC_VERSION: $petsc_version"
echo "PETSC_VERSION=$petsc_version" >> $GITHUB_ENV
- name: Set GCC-12 as default
run: |
echo "CC=gcc-12" >> $GITHUB_ENV
echo "CXX=g++-12" >> $GITHUB_ENV
- name: Cache PETSc Installation
id: cache-petsc
uses: actions/cache@v3
uses: actions/cache@v4
with:
path: ${{ env.PETSC_CACHE_FOLDER }}
key: ubuntu-${{ env.UBUNTU_VERSION }}-petsc-${{ env.PETSC_VERSION }}-${{ env.COMPILER }}-${{ env.ARCHITECTURE }}
restore-keys: ubuntu-${{ env.UBUNTU_VERSION }}-petsc-${{ env.PETSC_VERSION }}-${{ env.COMPILER }}-${{ env.ARCHITECTURE }}

- name: Install Dependencies
run: |
sudo apt-get update
sudo apt-get -y install \
gcc-12 \
g++-12 \
cmake \
openmpi-bin \
openmpi-common \
libopenmpi-dev \
libmetis-dev \
libparmetis-dev \
libeigen3-dev \
python3 \
python3-numpy \
libopenblas-base \
libopenblas-dev \
liblua5.3-0 \
liblua5.3-dev
key: ubuntu-${{ env.UBUNTU_VERSION }}-petsc-${{ env.PETSC_VERSION }}-${{ matrix.compiler.cc }}-${{ env.ARCHITECTURE }}
restore-keys: ubuntu-${{ env.UBUNTU_VERSION }}-petsc-${{ env.PETSC_VERSION }}-${{ matrix.compiler.cc }}-${{ env.ARCHITECTURE }}

- name: Install PETSc
if: steps.cache-petsc.outputs.cache-hit != 'true'
Expand All @@ -72,8 +70,8 @@ jobs:
cd petsc-${version}
./configure --with-fortran-bindings=0 --with-debugging=0 \
--with-memalign=32 --with-64-bit-indices \
CC=mpicc CXX=mpicxx FC=mpif90 --prefix=${install_dir}
COPTFLAGS="-g -O3" CXXOPTFLAGS="-g -O3"
--with-cc=${{ matrix.compiler.cc }} --with-cxx=${{ matrix.compiler.cxx }} --with-fc=0 --prefix=${install_dir} \
--COPTFLAGS="-g -O3" --CXXOPTFLAGS="-g -O3" --with-mpi-dir=/usr/lib/x86_64-linux-gnu/openmpi
make PETSC_DIR=`pwd` PETSC_ARCH=${{ env.ARCHITECTURE }} -j$(nproc)
make PETSC_DIR=`pwd` PETSC_ARCH=${{ env.ARCHITECTURE }} install
cd ..
Expand All @@ -82,7 +80,12 @@ jobs:
run: |
set -euo pipefail
mkdir build && cd build
cmake .. -DCMAKE_PREFIX_PATH=${{ env.PETSC_CACHE_FOLDER }} -DDOMAIN_DIMENSION=3 -DPOLYNOMIAL_DEGREE=2
cmake .. \
-DCMAKE_PREFIX_PATH=${{ env.PETSC_CACHE_FOLDER }} \
-DDOMAIN_DIMENSION=3 \
-DPOLYNOMIAL_DEGREE=2 \
-DCMAKE_C_COMPILER=${{ matrix.compiler.cc }} \
-DCMAKE_CXX_COMPILER=${{ matrix.compiler.cxx }}
make -j$(nproc)
make test
16 changes: 8 additions & 8 deletions app/form/SeasQDDiscreteGreenOperator.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -245,7 +245,7 @@ PetscInt SeasQDDiscreteGreenOperator::create_discrete_greens_function() {
CHKERRTHROW(MatCreateDense(comm, m, n, PETSC_DECIDE, PETSC_DECIDE, nullptr, &G_));
CHKERRTHROW(MatSetBlockSizes(G_, m_bs, n_bs));
CHKERRTHROW(MatGetSize(G_, &M, &N));
CHKERRTHROW(PetscPrintf(comm, "Green's function operator size: %D x %D\n", M, N));
CHKERRTHROW(PetscPrintf(comm, "Green's function operator size: %ld x %ld\n", M, N));

S_ = std::make_unique<PetscVector>(slip_block_size, num_local_elements, comm);
t_boundary_ = std::make_unique<PetscVector>(m_bs, num_local_elements, comm);
Expand Down Expand Up @@ -281,7 +281,7 @@ void SeasQDDiscreteGreenOperator::write_discrete_greens_operator(
"write_discrete_greens_operator():matrix %1.2e (sec)\n",
(double)(t1 - t0)));
CHKERRTHROW(PetscPrintf(PetscObjectComm((PetscObject)G_),
" status: computed %D / pending %D\n", current_gf, n_gf - current_gf));
" status: computed %ld / pending %ld\n", current_gf, n_gf - current_gf));

write_facet_labels_IS(mesh);
}
Expand Down Expand Up @@ -497,13 +497,13 @@ PetscInt SeasQDDiscreteGreenOperator::load_discrete_greens_operator(
MPI_Comm_size(comm, &commsize);
if ((PetscInt)commsize != commsize_checkpoint) {
CHKERRTHROW(PetscPrintf(comm,
"GF loaded was created with commsize = %d. Current commsize = %d. "
"GF loaded was created with commsize = %ld. Current commsize = %ld. "
"Repartitioning required.\n",
(int)commsize_checkpoint, (int)commsize));
(PetscInt)commsize_checkpoint, (PetscInt)commsize));
repartition_gfs_ = true;
} else {
CHKERRTHROW(PetscPrintf(
comm, "GF loaded was created with commsize matching current (%d).\n", (int)commsize));
comm, "GF loaded was created with commsize matching current (%ld).\n", (PetscInt)commsize));
repartition_gfs_ = false;
}

Expand All @@ -514,7 +514,7 @@ PetscInt SeasQDDiscreteGreenOperator::load_discrete_greens_operator(
CHKERRTHROW(PetscTime(&t1));
CHKERRTHROW(PetscPrintf(PetscObjectComm((PetscObject)G_),
"load_discrete_greens_operator() %1.2e (sec)\n", (double)(t1 - t0)));
CHKERRTHROW(PetscPrintf(PetscObjectComm((PetscObject)G_), " status: loaded %D / pending %D\n",
CHKERRTHROW(PetscPrintf(PetscObjectComm((PetscObject)G_), " status: loaded %ld / pending %ld\n",
current_gf, n_gf - current_gf));

if (repartition_gfs_) {
Expand Down Expand Up @@ -565,12 +565,12 @@ void SeasQDDiscreteGreenOperator::partial_assemble_discrete_greens_function(
auto ghost = scatter.template recv_prototype<double>(slip_block_size, ALIGNMENT);

CHKERRTHROW(PetscPrintf(PetscObjectComm((PetscObject)G_),
"partial_assemble_discrete_greens_function() [%D , %D)\n", start, N));
"partial_assemble_discrete_greens_function() [%ld , %ld)\n", start, N));
solve_time = 0.0;
for (PetscInt i = start; i < N; ++i) {

CHKERRTHROW(PetscPrintf(PetscObjectComm((PetscObject)G_),
"Computing Green's function %D / %D\n", i, N));
"Computing Green's function %ld / %ld\n", i, N));
sw.start();
CHKERRTHROW(VecZeroEntries(S_->vec()));
if (i >= nb_offset && i < nb_offset + m) {
Expand Down
14 changes: 7 additions & 7 deletions app/pc/eigdeflate.c
Original file line number Diff line number Diff line change
Expand Up @@ -173,7 +173,7 @@ PetscErrorCode PCDestroy_eigdeflate(PC pc) {
PetscFunctionReturn(0);
}

PetscErrorCode PCSetFromOptions_eigdeflate(PetscOptionItems* PetscOptionsObject, PC pc) {
PetscErrorCode PCSetFromOptions_eigdeflate(PC pc, PetscOptionItems* PetscOptionsObject) {
PC_eigdeflate* ctx;
PetscInt M;
PetscBool flg;
Expand Down Expand Up @@ -203,14 +203,14 @@ PetscErrorCode PCView_eigdeflate(PC pc, PetscViewer viewer) {
PC_eigdeflate* ctx = (PC_eigdeflate*)pc->data;

PetscFunctionBegin;
PetscViewerASCIIPrintf(viewer, "num. eigenvectors: %D\n", ctx->nev);
PetscViewerASCIIPrintf(viewer, "num. eigenvectors: %ld\n", ctx->nev);
PetscViewerASCIIPrintf(viewer, "emin: %+1.4e\n", ctx->e_min);
PetscViewerASCIIPrintf(viewer, "emax: %+1.4e\n", ctx->e_max);
PetscViewerASCIIPrintf(viewer, "Randomized eigenvalue calculation\n", ctx->nev_oversample);
PetscViewerASCIIPrintf(viewer, "over sampling: %D\n", ctx->nev_oversample);
PetscViewerASCIIPrintf(viewer, "power iterations: %D\n", ctx->power_its);
PetscViewerASCIIPrintf(viewer, "pre-smooth iterations: %D\n", ctx->npre);
PetscViewerASCIIPrintf(viewer, "post-smooth iterations: %D\n", ctx->npost);
PetscViewerASCIIPrintf(viewer, "Randomized eigenvalue calculation\n");
PetscViewerASCIIPrintf(viewer, "over sampling: %ld\n", ctx->nev_oversample);
PetscViewerASCIIPrintf(viewer, "power iterations: %ld\n", ctx->power_its);
PetscViewerASCIIPrintf(viewer, "pre-smooth iterations: %ld\n", ctx->npre);
PetscViewerASCIIPrintf(viewer, "post-smooth iterations: %ld\n", ctx->npost);
if (ctx->npre > 0 || ctx->npost > 0) {
PetscViewerASCIIPrintf(viewer, "optimal relaxation: %+1.4e\n", ctx->alpha);
}
Expand Down
8 changes: 4 additions & 4 deletions app/pc/reig_aux.c
Original file line number Diff line number Diff line change
Expand Up @@ -1093,8 +1093,8 @@ PetscErrorCode RandEigsMax(Mat A,PetscInt k,PetscInt o,PetscInt pits,PetscRandom
comm = PetscObjectComm((PetscObject)A);
ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
if (M != N) SETERRQ2(comm,PETSC_ERR_SUP,"Only valid for square matrices, found M = %D, N = %D\n",M,N);
if (k+o > M) SETERRQ4(comm,PETSC_ERR_SUP,"Random matrix has %D + %D = %D columns. Max num. columns is %D\n",k,o,k+o,M);
if (M != N) SETERRQ2(comm,PETSC_ERR_SUP,"Only valid for square matrices, found M = %ld, N = %ld\n",M,N);
if (k+o > M) SETERRQ4(comm,PETSC_ERR_SUP,"Random matrix has %ld + %ld = %ld columns. Max num. columns is %ld\n",k,o,k+o,M);
ierr = RandEigsMax_3_InPlace(A,k,o,pits,prand,_eigs,_V);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
Expand All @@ -1112,8 +1112,8 @@ PetscErrorCode RandEigsMin(KSP ksp,PetscInt k,PetscInt o,PetscInt pits,PetscRand
comm = PetscObjectComm((PetscObject)A);
ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
if (M != N) SETERRQ2(comm,PETSC_ERR_SUP,"Only valid for square matrices, found M = %D, N = %D\n",M,N);
if (k+o > M) SETERRQ4(comm,PETSC_ERR_SUP,"Random matrix has %D + %D = %D columns. Max num. columns is %D\n",k,o,k+o,M);
if (M != N) SETERRQ2(comm,PETSC_ERR_SUP,"Only valid for square matrices, found M = %ld, N = %ld\n",M,N);
if (k+o > M) SETERRQ4(comm,PETSC_ERR_SUP,"Random matrix has %ld + %ld = %ld columns. Max num. columns is %ld\n",k,o,k+o,M);
ierr = RandEigsMin_3_InPlace(ksp,k,o,pits,prand,_eigs,_V);CHKERRQ(ierr);

PetscFunctionReturn(0);
Expand Down
1 change: 1 addition & 0 deletions src/util/Schema.h
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,7 @@
#include <stdexcept>
#include <toml.hpp>

#include <algorithm>
#include <charconv>
#include <functional>
#include <limits>
Expand Down

0 comments on commit 82d2ee5

Please sign in to comment.