[petsc-users] one compilation error in PETSc-dev with enabling GPU and complex number

recrusader recrusader at gmail.com
Tue Feb 7 10:19:20 CST 2012


Dear Matt,

I got the help from Thrust developers.
The problem is
in PETSc, we define CUSPARRAY using cusp::array1d<PetscScalar,
cusp::device_memory>. Generaly, PetscScalar is std::complex.
However, the CUSPARRAY array cannot be decorated at GPU side. We need
cusp::array1d<cusp::complex, cusp::device_memory>.
Is there any simple method to process this in PETSc?

Thank you very much.

Best,
Yujie

On Sun, Jan 29, 2012 at 1:20 PM, Matthew Knepley <knepley at gmail.com> wrote:

> On Sun, Jan 29, 2012 at 1:05 PM, recrusader <recrusader at gmail.com> wrote:
>
>> Thank you very much, Matt,
>>
>> You mean the headers of the simple codes, I further simply the codes as
>
>
> This is a question for the CUSP mailing list.
>
>   Thanks,
>
>      Matt
>
>
>> "
>> #include <cusp/blas.h>
>> #include <cusp/array1d.h>
>>
>> int main(void)
>> {
>>    cusp::array1d<std::complex<double>, cusp::host_memory> *x;
>>
>>    x=new cusp::array1d<std::complex<double>, cusp::host_memory>(2,0.0);
>>
>>    std::complex<double> alpha(1,2.0);
>>    cusp::blas::scal(*x,alpha);
>>
>>    return 0;
>> }"
>>
>> I got the same compilation results "
>> login1$ nvcc gputest.cu -o gputest
>> /opt/apps/cuda/4.0/cuda/bin/../include/cusp/detail/blas.inl(134):
>> warning: calling a __host__ function from a __host__ __device__
>> function is not allowed
>>          detected during:
>>            instantiation of "void
>> cusp::blas::detail::SCAL<T>::operator()(T2 &) [with
>> T=std::complex<double>, T2=std::complex<double>]"
>> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/host/for_each.inl(37):
>> here
>>            instantiation of "InputIterator
>> thrust::detail::host::for_each(InputIterator, InputIterator,
>> UnaryFunction) [with
>> InputIterator=thrust::detail::normal_iterator<std::complex<double> *>,
>> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>>
>> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/dispatch/for_each.h(46):
>> here
>>            instantiation of "InputIterator
>> thrust::detail::dispatch::for_each(InputIterator, InputIterator,
>> UnaryFunction, thrust::host_space_tag) [with
>> InputIterator=thrust::detail::normal_iterator<std::complex<double> *>,
>> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/for_each.inl(51):
>> here
>>            instantiation of "InputIterator
>> thrust::detail::for_each(InputIterator, InputIterator, UnaryFunction)
>> [with InputIterator=thrust::detail::normal_iterator<std::complex<double>
>> *>, UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/for_each.inl(67):
>> here
>>            instantiation of "void thrust::for_each(InputIterator,
>> InputIterator, UnaryFunction) [with
>> InputIterator=thrust::detail::normal_iterator<std::complex<double> *>,
>> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> (367): here
>>            instantiation of "void
>> cusp::blas::detail::scal(ForwardIterator, ForwardIterator, ScalarType)
>> [with ForwardIterator=thrust::detail::normal_iterator<std::complex<double>
>> *>, ScalarType=std::complex<double>]"
>> (748): here
>>            instantiation of "void cusp::blas::scal(Array &,
>> ScalarType) [with Array=cusp::array1d<std::complex<double>,
>> cusp::host_memory>, ScalarType=std::complex<double>]"
>> gputest.cu(25): here
>>
>> "
>>
>> Thanks a lot.
>>
>> Best,
>> Yujie
>>
>> On 1/29/12, Matthew Knepley <knepley at gmail.com> wrote:
>> > On Sun, Jan 29, 2012 at 12:53 PM, recrusader <recrusader at gmail.com>
>> wrote:
>> >
>> >> Dear PETSc developers,
>> >>
>> >> With your help, I can successfully PETSc-deve with enabling GPU and
>> >> complex number.
>> >> However, when I compiled the codes, I met some errors. I also tried to
>> >> use simple codes to realize the same function. However, the errors
>> >> disappear. One example is as follows:
>> >>
>> >> for the function "VecScale_SeqCUSP"
>> >> "#undef __FUNCT__
>> >> #define __FUNCT__ "VecScale_SeqCUSP"
>> >> PetscErrorCode VecScale_SeqCUSP(Vec xin, PetscScalar alpha)
>> >> {
>> >>  CUSPARRAY      *xarray;
>> >>  PetscErrorCode ierr;
>> >>
>> >>  PetscFunctionBegin;
>> >>  if (alpha == 0.0) {
>> >>    ierr = VecSet_SeqCUSP(xin,alpha);CHKERRQ(ierr);
>> >>  } else if (alpha != 1.0) {
>> >>    ierr = VecCUSPGetArrayReadWrite(xin,&xarray);CHKERRQ(ierr);
>> >>    try {
>> >>      cusp::blas::scal(*xarray,alpha);
>> >>    } catch(char* ex) {
>> >>      SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
>> >>    }
>> >>    ierr = VecCUSPRestoreArrayReadWrite(xin,&xarray);CHKERRQ(ierr);
>> >>  }
>> >>  ierr = WaitForGPU();CHKERRCUSP(ierr);
>> >>  ierr = PetscLogFlops(xin->map->n);CHKERRQ(ierr);
>> >>  PetscFunctionReturn(0);
>> >> } "
>> >>
>> >> When I compiled PETSc-dev, I met the following errors:
>> >> " /opt/apps/cuda/4.0/cuda/include/cusp/detail/blas.inl(134): warning:
>> >> calling a __host__ function from a __host__ __device__ function is not
>> >> allowed
>> >>          detected during:
>> >>            instantiation of "void
>> >> cusp::blas::detail::SCAL<T>::operator()(T2 &) [with
>> >> T=std::complex<double>, T2=PetscScalar]"
>> >>
>> /opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/for_each.inl(72):
>> >> here
>> >>            instantiation of "void
>> >> thrust::detail::device::cuda::for_each_n_closure<RandomAccessIterator,
>> >> Size, UnaryFunction>::operator()() [with
>> >>
>> >>
>> RandomAccessIterator=thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> Size=long,
>> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >>
>> >>
>> /opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/detail/launch_closure.inl(51):
>> >> here
>> >>            instantiation of "void
>> >>
>> >>
>> thrust::detail::device::cuda::detail::launch_closure_by_value(NullaryFunction)
>> >> [with
>> >>
>> NullaryFunction=thrust::detail::device::cuda::for_each_n_closure<thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> long, cusp::blas::detail::SCAL<std::complex<double>>>]"
>> >>
>> >>
>> /opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/detail/launch_closure.inl(71):
>> >> here
>> >>            instantiation of "size_t
>> >>
>> >>
>> thrust::detail::device::cuda::detail::closure_launcher_base<NullaryFunction,
>> >> launch_by_value>::block_size_with_maximal_occupancy(size_t) [with
>> >>
>> >>
>> NullaryFunction=thrust::detail::device::cuda::for_each_n_closure<thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> long, cusp::blas::detail::SCAL<std::complex<double>>>,
>> >> launch_by_value=true]"
>> >>
>> >>
>> /opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/detail/launch_closure.inl(136):
>> >> here
>> >>            instantiation of "thrust::pair<size_t, size_t>
>> >>
>> >>
>> thrust::detail::device::cuda::detail::closure_launcher<NullaryFunction>::configuration_with_maximal_occupancy(Size)
>> >> [with
>> >>
>> NullaryFunction=thrust::detail::device::cuda::for_each_n_closure<thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> long, cusp::blas::detail::SCAL<std::complex<double>>>, Size=long]"
>> >>
>> >>
>> /opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/detail/launch_closure.inl(145):
>> >> here
>> >>            [ 6 instantiation contexts not shown ]
>> >>            instantiation of "InputIterator
>> >> thrust::detail::dispatch::for_each(InputIterator, InputIterator,
>> >> UnaryFunction, thrust::device_space_tag) [with
>> >>
>> >>
>> InputIterator=thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >> /opt/apps/cuda/4.0/cuda/include/thrust/detail/for_each.inl(51): here
>> >>            instantiation of "InputIterator
>> >> thrust::detail::for_each(InputIterator, InputIterator, UnaryFunction)
>> >> [with
>> >>
>> InputIterator=thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >> /opt/apps/cuda/4.0/cuda/include/thrust/detail/for_each.inl(67): here
>> >>            instantiation of "void thrust::for_each(InputIterator,
>> >> InputIterator, UnaryFunction) [with
>> >>
>> >>
>> InputIterator=thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >> (367): here
>> >>            instantiation of "void
>> >> cusp::blas::detail::scal(ForwardIterator, ForwardIterator, ScalarType)
>> >> [with
>> >>
>> ForwardIterator=thrust::detail::normal_iterator<thrust::device_ptr<PetscScalar>>,
>> >> ScalarType=std::complex<double>]"
>> >> (748): here
>> >>            instantiation of "void cusp::blas::scal(Array &,
>> >> ScalarType) [with Array=cusp::array1d<PetscScalar,
>> >> cusp::device_memory>, ScalarType=std::complex<double>]"
>> >> veccusp.cu(1185): here
>> >>
>> >>
>> >>
>> /opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/detail/launch_closure.inl(51):
>> >> error: a value of type "int" cannot be assigned to an entity of type
>> >> "_ZNSt7complexIdE9_ComplexTE"
>> >>
>> >> "
>> >> However, I further realize simiar codes as
>> >> "
>> >> #include <thrust/version.h>
>> >> #include <cusp/version.h>
>> >> #include <iostream>
>> >> #include <cusp/blas.h>
>> >> #include <cusp/array1d.h>
>> >> #include <complex>
>> >>
>> >> int main(void)
>> >> {
>> >>    cusp::array1d<std::complex<double>, cusp::host_memory> *x;
>> >>
>> >>    x=new cusp::array1d<std::complex<double>, cusp::host_memory>(2,0.0);
>> >>
>> >>    std::complex<double> alpha(1,2.0);
>> >>    cusp::blas::scal(*x,alpha);
>> >>
>> >>    return 0;
>> >> }
>> >> "
>> >>
>> >> When I complied it using  "nvcc gputest.cu -o gputest", I only meet
>> >> warning information as follows:
>> >> "
>> >> /opt/apps/cuda/4.0/cuda/bin/../include/cusp/detail/blas.inl(134):
>> >> warning: calling a __host__ function from a __host__ __device__
>> >> function is not allowed
>> >>          detected during:
>> >>            instantiation of "void
>> >> cusp::blas::detail::SCAL<T>::operator()(T2 &) [with
>> >> T=std::complex<double>, T2=std::complex<double>]"
>> >>
>> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/host/for_each.inl(37):
>> >> here
>> >>            instantiation of "InputIterator
>> >> thrust::detail::host::for_each(InputIterator, InputIterator,
>> >> UnaryFunction) [with
>> >> InputIterator=thrust::detail::normal_iterator<std::complex<double> *>,
>> >> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >>
>> >>
>> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/dispatch/for_each.h(46):
>> >> here
>> >>            instantiation of "InputIterator
>> >> thrust::detail::dispatch::for_each(InputIterator, InputIterator,
>> >> UnaryFunction, thrust::host_space_tag) [with
>> >> InputIterator=thrust::detail::normal_iterator<std::complex<double> *>,
>> >> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/for_each.inl(51):
>> >> here
>> >>            instantiation of "InputIterator
>> >> thrust::detail::for_each(InputIterator, InputIterator, UnaryFunction)
>> >> [with
>> InputIterator=thrust::detail::normal_iterator<std::complex<double>
>> >> *>, UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >> /opt/apps/cuda/4.0/cuda/bin/../include/thrust/detail/for_each.inl(67):
>> >> here
>> >>            instantiation of "void thrust::for_each(InputIterator,
>> >> InputIterator, UnaryFunction) [with
>> >> InputIterator=thrust::detail::normal_iterator<std::complex<double> *>,
>> >> UnaryFunction=cusp::blas::detail::SCAL<std::complex<double>>]"
>> >> (367): here
>> >>            instantiation of "void
>> >> cusp::blas::detail::scal(ForwardIterator, ForwardIterator, ScalarType)
>> >> [with
>> ForwardIterator=thrust::detail::normal_iterator<std::complex<double>
>> >> *>, ScalarType=std::complex<double>]"
>> >> (748): here
>> >>            instantiation of "void cusp::blas::scal(Array &,
>> >> ScalarType) [with Array=cusp::array1d<std::complex<double>,
>> >> cusp::host_memory>, ScalarType=std::complex<double>]"
>> >> gputest.cu(25): here
>> >>
>> >> "
>> >> There are not errors like
>> >>
>> >>
>> "/opt/apps/cuda/4.0/cuda/include/thrust/detail/device/cuda/detail/launch_closure.inl(51):
>> >> error: a value of type "int" cannot be assigned to an entity of type
>> >> "_ZNSt7complexIdE9_ComplexTE" "
>> >>
>> >> Furthermore, the warning information is also different between
>> >> PETSc-dev and simple codes.
>> >>
>> >> Could you give me some suggestion for this errors? Thank you very much.
>> >>
>> >
>> > The headers are complicated to get right. The whole point of what we
>> did is
>> > to give a way to use GPU
>> > simply through the existing PETSc linear algebra interface.
>> >
>> >    Matt
>> >
>> >
>> >> Best,
>> >> Yujie
>> >>
>> >
>> >
>> >
>> > --
>> > 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
>> >
>>
>
>
>
> --
> 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
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120207/7279e2c2/attachment-0001.htm>


More information about the petsc-users mailing list