[petsc-users] compiling petsc.h and cusp's headers

Satish Balay balay at mcs.anl.gov
Tue Oct 4 21:32:47 CDT 2011


How about the following?

Satish

------------
balay at bb30:~/ex>cat makefile 

CFLAGS           =
FFLAGS           =
CPPFLAGS         =
FPPFLAGS         =
CLEANFILES       =

include ${PETSC_DIR}/conf/variables
include ${PETSC_DIR}/conf/rules

ex1: ex1.o  chkopts
     -${CLINKER} -o ex1 ex1.o  ${PETSC_KSP_LIB}
     ${RM} ex1.o

balay at bb30:~/ex>cat ex1.cu
#include "petsc.h"
#undef VecType

#include<iostream>
#include <cusp/krylov/cg.h>

int main(int argc,char **argv)
{
  PetscErrorCode ierr;
  ierr = PetscInitialize(&argc,&argv,(char *)0,(char*)0);CHKERRQ(ierr);
  std::cout<<"Hello World!"<<std::endl;
  ierr = PetscFinalize();
  return 0;
}
balay at bb30:~/ex>make ex1
nvcc -g -arch=sm_13  -c --compiler-options="-Wall -Wwrite-strings -Wno-strict-aliasing -Wno-unknown-pragmas -g3   -I/home/balay/petsc-dev/include -I/home/balay/petsc-dev/arch-cuda-double/include -I/usr/local/cuda/include    -D__INSDIR__="  ex1.cu
/home/balay/petsc-dev/arch-cuda-double/bin/mpicc -Wall -Wwrite-strings -Wno-strict-aliasing -Wno-unknown-pragmas -g3  -o ex1 ex1.o  -L/home/balay/petsc-dev/arch-cuda-double/lib  -lpetsc -lX11 -lpthread -Wl,-rpath,/usr/local/cuda/lib64 -L/usr/local/cuda/lib64 -lcufft -lcublas -lcudart -llapack -lblas -lm -lmpichcxx -lstdc++ -ldl 
/bin/rm -f ex1.o
balay at bb30:~/ex>./ex1
Hello World!
balay at bb30:~/ex>


On Tue, 4 Oct 2011, Shiyuan wrote:

> *
> *>*  Are you using the PETSc makefiles to compile it? Or did you make up that
> *>* compile line yourself? Why not just use the PETSc makefile to compile it?
> 
> 
> 
> *Actually, the compiling command I usedd is generated by
> PETSs's makefiles.* *But I am looking for something like "make
> getincludedirs; make getpetscflags;
> make getlinklibs" which can give me the compiling flags Petsc uses to
> feed nvcc. Are there anything
> like that?
> *
> 
> *>*
> *>*  Presumably cusp/thrust has some pattern for what include files need to be
> *>* included in what order. If you do not know the pattern then why not start by
> *>* copying what PETSc does, since that compiles? For example
> *>* cuspvecimpl.h has
> *>*
> *>* #include <cublas.h>
> *>* #include <cusp/blas.h>
> *>* #include <thrust/device_vector.h>
> *>* #include <thrust/iterator/constant_iterator.h>
> *>* #include <thrust/transform.h>
> *>* #include <thrust/iterator/permutation_iterator.h>
> *>*
> *>* while cuspmatimpl.h has
> *>*
> *>* #include <../src/vec/vec/impls/seq/seqcusp/cuspvecimpl.h>
> *>* #include <cusp/csr_matrix.h>
> *>* #include <cusp/multiply.h>
> *>* /*for MatCreateSeqAIJCUSPFromTriple*/
> *>* #include <cusp/coo_matrix.h>
> *>*
> *>* so start with at least those.   Generally just picking a couple of random
> *>* include files from some C++ package (like cusp/thrust)  and including them
> *>* won't work.
> 
> 
> I did check the documentation of cusp, and I don't see a requirement about
> **the order for including headers. It doesn't sound reasonable to
> require programmers to enforce the order when the programmers
> uses only functions from one library. Does it?( I think there is an
> issue I don't see, please excuse me for my ignorance). And I did
> try a few experiment where I picked some cusp examples and changed the
> orders of including headers, and they do work.
> 
> 
> 
> 
> *>I believe that your problem is the clash of VecType between CUSP and PETSc.
> >In aijcusp.cu, we have
> 
> >#undef VecType
> >#include "../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h"
> 
>  >   Matt
> 
> That explains it. But I still don't get the solution. If I want to
> define a function where I call some functions provided by Petsc and
> some functions  provided by cusp and thrust,
> In what order should I include the headers? For example, if I want
> Mat, Vec and cusp::krylov::cg.
>  If I look at the aijcusp.cu where MatMult_SeqAIJCusp is defined, it
> has a long list of headers,
> 
> 
> #include "petscconf.h"
> PETSC_CUDA_EXTERN_C_BEGIN
> #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
> #include "petscbt.h"
> #include "../src/vec/vec/impls/dvecimpl.h"
> #include "private/vecimpl.h"
> PETSC_CUDA_EXTERN_C_END
> #undef VecType
> #include "../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h"
> 
> 
> #ifdef PETSC_HAVE_TXPETSCGPU
> 
> #include "csr_matrix_data.h"
> #include "csr_matrix_data_gpu.h"
> #include "csr_tri_solve_gpu.h"
> #include "csr_tri_solve_gpu_level_scheduler.h"
> #include "csr_spmv_inode.h"
> #include <algorithm>
> #include <vector>
> #include <string>
> #include <thrust/sort.h>
> #include <thrust/fill.h>
> 
> 
> But looking at this list of headers, I cannot see what files do I
> really need and what else is missing if I only want Mat, Vec and
> cusp::krylov::cg.
> Is there a less painful way? Thanks.
> 
> 
> Shiyuan
> ****
> 
> ---------------------------------------------------------------------------------------------------------------------------------------------------------------------
> **>* On Oct 4, 2011, at 2:07 PM, Shiyuan wrote:
> *>*
> *>* > Hi,
> *>* >    I cannot compiling a source file with petsc.h and cusp headers like
> *>* the following:
> *>* >
> *>* > #include"petsc.h"
> *>* > #include<iostream>
> *>* > #include "cusp/csr_matrix.h"
> *>* > #include <cusp/krylov/cg.h>
> *>* > int main(void){
> *>* >     std::cout<<"Hello World!"<<std::endl;
> *>* >     return 0;
> *>* > }
> *>* >
> *>* > I use the following  to compile it:
> *>* >
> *>* > nvcc  -m64 -O -arch=sm_13  -c --compiler-options="-Wall -Wwrite-strings
> *>* -Wno-strict-aliasing -Wno-unknown-pragmas -O
> *>* -I/home/guest/sgu1/softwares/petsc-dev/include
> *>* -I/home/guest/sgu1/softwares/petsc-dev/helena-cxx-nompi-os64-release/include
> *>* -I/usr/local/cuda/include -I/home/guest/sgu1/softwares/cusp-library/
> *>* -I/home/guest/sgu1/softwares/thrust/
> *>* -I/home/guest/sgu1/softwares/petsc-dev/include/mpiuni  -DCUDA=1
> *>* -I/home/guest/sgu1/softwares/slepc-dev
> *>* -I/home/guest/sgu1/softwares/slepc-dev/helena-cxx-nompi-os64-release/include
> *>* -I/home/guest/sgu1/softwares/slepc-dev/include
> *>* -I/home/guest/sgu1/softwares/ImageMagick-6.7.2/include/ImageMagick
> *>*  -D__INSDIR__= -I/home/guest/sgu1/softwares/slepc-dev
> *>* -I/home/guest/sgu1/softwares/slepc-dev/helena-cxx-nompi-os64-release/include
> *>* -I/home/guest/sgu1/softwares/slepc-dev/include" MyKspTmp4.cu -o MyKspTmp4.o
> *>* >
> *>* > and it gives me error:
> *>* >
> *>* /usr/local/cuda/include/thrust/detail/device/cuda/detail/b40c/vector_types.h(37):
> *>* error: expected an identifier
> *>* >
> *>* > I think I did not do it in a correct way. I took a a look at some of
> *>* petsc's .cu files, they usually include a  long list of headers. Can petsc.h
> *>* directly included in a source code with cusp headers? What's the correct way
> *>* to do it?  Does the makefile system in petsc provide the command or flags to
> *>* do that? How can I obtain those flags?
> *>* >
> *>* > And after I compile it and obtain .o files, can I use g++ to link the .o
> *>* files with petcs libs and other .o files as usual, i.e.,
> *>* >     g++ -o $(BIN_DIR)/$@ $(CPPFLAGS) $@.s $(objects) $(PETSC_LIB)
> *>* $(OTHERS)
> *>* > or I need to do some extra steps and some specific linking flags to do
> *>* that?
> *>* > Thanks.
> *>* >
> *>* > Shiyuan
> *>* >
> *>* >
> *>* >
> *>*
> *>*
> *
> 



More information about the petsc-users mailing list