[petsc-dev] petsc/master: unable to link in C++ with last night PETSC_WITH_EXTERNAL_LIB variable changes
Smith, Barry F.
bsmith at mcs.anl.gov
Sat Feb 10 12:02:45 CST 2018
Eric,
How difficult would it be for you guys to switch to doing things "correctly" and not mixing up C++ compilers/linkers?
Barry
> On Feb 10, 2018, at 11:51 AM, Éric Chamberland <Eric.Chamberland at giref.ulaval.ca> wrote:
>
> Ok guys!
>
> You found it (again...). We are no using mpic++ to link!!!
>
> We are not always using the same compiler that we used to compile petsc+mpi, so we do not use mpic++ to compile (excepted for intelmpi).
>
> For example, I compile MPI+Petsc with g++, but I am working with clang++ and may switch to g++ for my day to day work. This is perfectly working since... long long time ago (I won't reveal my age... ;) )
>
> I understand you choice, but your PETSC_WITH_EXTERNAL_LIB variable was helping me to painlessly find the good linking options.
>
> I can work to manually add it when it is relevant into our compilation flags...
>
> How did you managed to extract the libs from all the different flavors of MPI (mpich, openmpi, intelmpi) so I can now reproduce what you did before?
>
> (maybe this should have been my question... )
>
> Thanks again,
>
> Eric
>
> Le 18-02-10 à 12:37, Smith, Barry F. a écrit :
>> Eric,
>>
>> We need the entire link line.
>>
>> What linker are you using C, or C++? This is important.
>>
>> Do you have dependencies on MPI C++ symbols? In other words, are you using C++ MPI bindings?
>>
>> I cannot explain why all the libraries you listed would disappear from PETSC_WITH_EXTERNAL_LIB but we did recently make a change to manual pass under some circumstances less system (including MPI libraries) explicitly since they are already usually passed by the linker. This may be causing your difficulties.
>>
>> Eagerly awaiting your reply.
>>
>> Barry
>>
>>
>>> On Feb 10, 2018, at 11:29 AM, Éric Chamberland <Eric.Chamberland at giref.ulaval.ca> wrote:
>>>
>>> Hi Matthew,
>>>
>>> Yes, I heard that the MPI C++ API has been deprecated.
>>>
>>> Yes, the mpi_cxx is now missing. Our link line is formed mainly with PETSC_WITH_EXTERNAL_LIB variable that is now:
>>>
>>> PETSC_WITH_EXTERNAL_LIB = -L/opt/petsc-master_debug/lib -Wl,-rpath,/opt/petsc-master_debug/lib -L/opt/petsc-master_debug/lib -L/opt/intel/composer_xe_2015.2.164/mkl/lib/intel64 -Wl,-rpath,/opt/intel/composer_xe_2015.2.164/mkl/lib/intel64 -Wl,-rpath,/opt/openmpi-1.10.2/lib -L/opt/openmpi-1.10.2/lib -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.8.5 -L/usr/lib/gcc/x86_64-redhat-linux/4.8.5 -Wl,-rpath,/opt/intel/composer_xe_2015.2.164/compiler/lib/intel64 -L/opt/intel/composer_xe_2015.2.164/compiler/lib/intel64 -Wl,-rpath,/opt/intel/composer_xe_2015.2.164/ipp/lib/intel64 -L/opt/intel/composer_xe_2015.2.164/ipp/lib/intel64 -Wl,-rpath,/opt/intel/composer_xe_2015.2.164/tbb/lib/intel64/gcc4.4 -L/opt/intel/composer_xe_2015.2.164/tbb/lib/intel64/gcc4.4 -lpetsc -lsuperlu -lsuperlu_dist -lHYPRE -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lmkl_scalapack_lp64 -lmkl_blacs_openmpi_lp64 -lml -lumfpack -lklu -lcholmod -lbtf -lccolamd -lcolamd -lcamd -lamd -lsuitesparseconfig -lmkl_intel_lp64 -lmkl_core -lmkl_intel_thread -lmkl_blacs_intelmpi_lp64 -liomp5 -ldl -lpthread -lparmetis -lmetis -lptesmumps -lptscotch -lptscotcherr -lesmumps -lscotch -lscotcherr -lm -lX11 -lstdc++ -ldl -lmpi_usempi -lmpi_mpifh -lmpi -lgfortran -lm -lgfortran -lm -lgcc_s -lquadmath -lpthread -lrt -lm -lpthread -lz -lstdc++ -ldl
>>> We add some stuff to this (our own compiled libs), but nothing related to MPI or PETSc since we used to rely on PETSC_WITH_EXTERNAL_LIB for all our diffrement environments (think about different petsc version, MPI libs and compilers: clang, icc, g++) and it used to work until yesterday changes...
>>>
>>> As you can see in the diff, this "block" of libraries have been removed from PETSC_WITH_EXTERNAL_LIB: (the -- is part of the diff output):
>>> --ldl
>>> --lmpi_cxx
>>> --lmpi
>>> --lstdc++
>>> --lm
>>> --lgcc_s
>>> --lpthread
>>>
>>> and into this block there was the "-lmpi_cxx" that we need...
>>>
>>> I could send you our whole line of link, but the error is into this small change introduced yesterday into master...
>>>
>>> Thanks a lot!
>>>
>>> Eric
>>>
>>> Le 18-02-10 à 10:34, Matthew Knepley a écrit :
>>>> On Sat, Feb 10, 2018 at 9:42 AM, Éric Chamberland <Eric.Chamberland at giref.ulaval.ca> wrote:
>>>> Hi,
>>>>
>>>> we used to link our c++ code with PETSc using PETSC_WITH_EXTERNAL_LIB variable defined in $PETSC_DIR/lib/petsc/conf/petscvariables and everything was fine until this night.
>>>>
>>>> It seems some libs have vanished from this variable, see the diff here:
>>>>
>>>> -lptscotcherr
>>>> -lesmumps
>>>> -lscotch
>>>> -lscotcherr
>>>> -lm
>>>> -lX11
>>>> +-lstdc++
>>>> -ldl
>>>> -lmpi_usempi
>>>> -lmpi_mpifh
>>>> -lmpi
>>>> -lgfortran
>>>> -lm
>>>> -lgfortran
>>>> -lm
>>>> -lgcc_s
>>>> -lquadmath
>>>> -lpthread
>>>> --ldl
>>>> --lmpi_cxx
>>>> --lmpi
>>>> --lstdc++
>>>> --lm
>>>> --lgcc_s
>>>> --lpthread
>>>> -lrt
>>>> -lm
>>>> -lpthread
>>>> -lz
>>>> +-lstdc++
>>>> -ldl
>>>>
>>>>
>>>> causing these errors at link phase for us:
>>>>
>>>>
>>>> /pmi/cmpbib/compilation_BIB_gcc_redhat_petsc-master_debug/COMPILE_AUTO/GIREF/obj/dev/StatistiqueMemoire.o: In function `MPI::Op::Init(void (*)(void const*, void*, int, MPI::Datatype const&), bool)':
>>>>
>>>> /opt/openmpi-1.10.2/include/
>>>> openmpi/ompi/mpi/cxx/op_inln.
>>>> h:122: undefined reference to `ompi_mpi_cxx_op_intercept'
>>>>
>>>> /pmi/cmpbib/compilation_BIB_gcc_redhat_petsc-master_debug/COMPILE_AUTO/GIREF/obj/dev/StatistiqueMemoire.o: In function `MPI::Intracomm::Create_graph(int, int const*, int const*, bool) const':
>>>>
>>>> /opt/openmpi-1.10.2/include/
>>>> openmpi/ompi/mpi/cxx/
>>>> intracomm.h:25: undefined reference to `MPI::Comm::Comm()'
>>>>
>>>> /pmi/cmpbib/compilation_BIB_gcc_redhat_petsc-master_debug/COMPILE_AUTO/GIREF/obj/dev/StatistiqueMemoire.o: In function `MPI::Intercomm::Merge(bool) const':
>>>>
>>>> /opt/openmpi-1.10.2/include/
>>>> openmpi/ompi/mpi/cxx/
>>>> intracomm_inln.h:23: undefined reference to `MPI::Comm::Comm()'
>>>>
>>>> /pmi/cmpbib/compilation_BIB_gcc_redhat_petsc-master_debug/COMPILE_AUTO/GIREF/obj/dev/StatistiqueMemoire.o: In function `MPI::Intracomm::Split(int, int) const':
>>>>
>>>> /opt/openmpi-1.10.2/include/
>>>> openmpi/ompi/mpi/cxx/
>>>> intracomm_inln.h:23: undefined reference to `MPI::Comm::Comm()'
>>>>
>>>> /pmi/cmpbib/compilation_BIB_gcc_redhat_petsc-master_debug/COMPILE_AUTO/GIREF/obj/dev/StatistiqueMemoire.o: In function `MPI::Intracomm::Create(MPI::Group const&) const':
>>>>
>>>> /opt/openmpi-1.10.2/include/
>>>> openmpi/ompi/mpi/cxx/
>>>> intracomm_inln.h:23: undefined reference to `MPI::Comm::Comm()'
>>>>
>>>> /pmi/cmpbib/compilation_BIB_gcc_redhat_petsc-master_debug/COMPILE_AUTO/GIREF/obj/dev/StatistiqueMemoire.o: In function `MPI::Intracomm::Clone() const':
>>>> Hi Eric,
>>>>
>>>> These symbols are all coming from -lmpi_cxx. I would note that I believe the MPI Forum has deprecated the C++ interface, so it will
>>>> eventually go away. However, lets fix this. In the configure log, I see that mpicxx -show has that library in it, so it seems that the link
>>>> is not being done with the C++ compiler. Can you send the whole link line?
>>>>
>>>> Thanks,
>>>>
>>>> Matt
>>>> Is this a normal and definitive change or an unwanted/unobserved bug?
>>>>
>>>> Thanks,
>>>>
>>>> Eric
>>>>
>>>> ps: here are the logs:
>>>>
>>>> this night:
>>>>
>>>> ---------
>>>>
>>>> http://www.giref.ulaval.ca/~cmpgiref/petsc-master-debug/2018.02.10.02h00m01s_configure.log
>>>> http://www.giref.ulaval.ca/~cmpgiref/petsc-master-debug/2018.02.10.02h00m01s_make.log
>>>>
>>>>
>>>>
>>>> a day before:
>>>> ------------
>>>>
>>>>
>>>> http://www.giref.ulaval.ca/~cmpgiref/petsc-master-debug/2018.02.09.02h00m02s_configure.log
>>>> http://www.giref.ulaval.ca/~cmpgiref/petsc-master-debug/2018.02.09.02h00m02s_make.log
>>>>
>>>>
>>>>
>>>>
>>>>
>>>>
>>>>
>>>> --
>>>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>>>> -- Norbert Wiener
>>>>
>>>> https://www.cse.buffalo.edu/~knepley/
>
>
More information about the petsc-dev
mailing list