[petsc-users] single precision vs double: strange behavior
Andreas Mang
andreas at ices.utexas.edu
Sat Mar 18 00:05:08 CDT 2017
Hey guys:
I was trying to run my code with single precision which resulted in strange behavior. The code essentially already breaks within the test case creation (computing some sinusoidal functions). I discovered that. e.g., norms and max and min values do not make sense. I created a simple test example that demonstrates strange behavior on my machine. I am trying to find the min and max values of a PETSc vector with numbers ranging from 0 to 15. I removed all non-essential compiler flags and dependencies.
If I compile PETSc with double precision I get
0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15
min 0 at 0
max 15 at 15
If I compile with single precision I get
0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15
min 15 at 8
max 14 at 7
The difference is that I add "--with-precision=single” to the configure command of PETSc. I use intel/17.0 with mpich2/3.1.4. The code can be found below. I also provide the configure options and the command line for compiling the code.
Thanks /Andreas
mpicxx -O3 -I./include -isystem/h1/andreas/code/develop/cold/external/libs/petsc_single/build/include -isystem/h1/andreas/code/develop/cold/external/libs/petsc_single/build/cxx_opt/include -c apps/checkcoldreg.cpp -o obj/checkcoldreg.o
Configure Options: --with-cc=mpicc --CFLAGS= --with-cxx=mpicxx --CXXFLAGS= --download-f2cblaslapack --with-debugging=1 --with-64-bit-indices --with-shared-libraries=0 --with-x=0 --with-fc=0 --with-precision=single
#include <iostream>
#include "petsc.h"
int main(int argc, char **argv) {
PetscErrorCode ierr;
PetscScalar maxval, minval;
PetscInt posmin, posmax;
Vec x; PetscScalar* p_x = NULL;
ierr = PetscInitialize(0, reinterpret_cast<char***>(NULL),
reinterpret_cast<char*>(NULL),
reinterpret_cast<char*>(NULL)); CHKERRQ(ierr);
PetscInt n = 16;
ierr = VecCreate(PETSC_COMM_WORLD, &x); CHKERRQ(ierr);
ierr = VecSetSizes(x, n, n); CHKERRQ(ierr);
ierr = VecSetFromOptions(x); CHKERRQ(ierr);
ierr = VecGetArray(x, &p_x); CHKERRQ(ierr);
for (PetscInt i = 0; i < n; ++i) {
p_x[i] = static_cast<PetscScalar>(i);
std::cout << p_x[i] << " ";
}
ierr = VecRestoreArray(x, &p_x); CHKERRQ(ierr);
std::cout << std::endl;
ierr = VecMin(x, &posmin, &minval); CHKERRQ(ierr);
ierr = VecMax(x, &posmax, &maxval); CHKERRQ(ierr);
std::cout << "min " << minval << " at " << posmin << std::endl;
std::cout << "max " << maxval << " at " << posmax << std::endl;
ierr = VecDestroy(&x); CHKERRQ(ierr);
ierr = PetscFinalize(); CHKERRQ(ierr);
return 0;
}
More information about the petsc-users
mailing list