[petsc-users] petscfe_type opencl, variable_coefficient field

David Fuentes fuentesdt at gmail.com
Mon Apr 21 22:30:50 CDT 2014


Hi,

Building off of snes/ex12,

I'm trying to use

-variable_coefficient field  with -petscfe_type opencl -mat_petscfe_type
 opencl

constant coefficients "-variable_coefficient none" works fine, but I am
getting NaN in my residual vector for spatially varying coefficients.

I tried changing N_t manually for N_comp_aux =1, also tried setting my
variable coefficient field to a constant at each quadrature point to see if
any difference.

Is there any other place I may look ?

https://bitbucket.org/petsc/petsc/src/fb3a4ffdea0449562fb0f58dad4b44552b6d3589/src/dm/dt/interface/dtfe.c?at=master#cl-4924

  if (useAux) {  ierr = PetscSNPrintfCount(string_tail, end_of_buffer
- string_tail,"    /* Load coefficients a_i for this cell */\n""    /*
TODO: This should not be N_t here, it should be N_bc*N_comp_aux */\n""
   a_i[tidx] = coefficientsAux[Goffset+batch*N_t+tidx];\n",
                &count);STRING_ERROR_CHECK("Message to short");  }



$ ./ireSolver -run_type full -dim 3 -petscspace_order 1
-mat_petscspace_order 0 -variable_coefficient field  -snes_type ksponly
 -snes_monitor -snes_converged_reason -ksp_converged_reason -ksp_rtol
1.e-12 -pc_type bjacobi -vtk geserviceStudy0025axt1.vtk -f meshIRElores.e
-sourcelandmark sourcelandmarks.vtk -targetlandmark targetlandmarks.vtk
-petscfe_type opencl -mat_petscfe_type  opencl -info -info_exclude
null,vec,mat -petscfe_num_blocks 2 -petscfe_num_batches 2
.
.
.
[0] PetscFEIntegrateResidual_OpenCL(): GPU layout grid(32,59,1)
block(8,1,1) with 2 batches
[0] PetscFEIntegrateResidual_OpenCL():  N_t: 8, N_cb: 2
[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[0]PETSC ERROR: Floating point exception
[0]PETSC ERROR: Vec entry at local location 120 is not-a-number or infinite
at end of function: Parameter number 3
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for
trouble shooting.
[0]PETSC ERROR: Petsc Development GIT revision: v3.4.4-4111-g7562b2e  GIT
Date: 2014-04-18 18:37:01 -0500
[0]PETSC ERROR: ./ireSolver on a arch-precise-gcc-4.6.3-dbg named maeda by
fuentes Mon Apr 21 22:23:25 2014
[0]PETSC ERROR: Configure options --with-shared-libraries
--with-clanguage=c++ --CFLAGS=-O0 --CXXFLAGS=-O0 --download-ctetgen
--download-triangle --with-debugging=yes
--with-blas-lapack-dir=/opt/apps/MKL/12.1/lib/intel64/
--with-exodusii-lib="[/usr/lib/x86_64-linux-gnu/libexoIIv2.so]"
--with-netcdf-dir=/usr/lib --with-hdf5-dir=/usr/lib --with-c2html=0
--with-exodusii-include=/usr/include --with-opencl-include=/usr/include/
--with-opencl-lib=/usr/lib/libOpenCL.so --download-viennacl=yes
[0]PETSC ERROR: #1 VecValidValues() line 34 in
/opt/apps/PETSc/petsc-usr/src/vec/vec/interface/rvector.c
[0]PETSC ERROR: #2 SNESComputeFunction() line 2096 in
/opt/apps/PETSc/petsc-usr/src/snes/interface/snes.c
[0]PETSC ERROR: #3 SNESSolve_KSPONLY() line 24 in
/opt/apps/PETSc/petsc-usr/src/snes/impls/ksponly/ksponly.c
[0]PETSC ERROR: #4 SNESSolve() line 3794 in
/opt/apps/PETSc/petsc-usr/src/snes/interface/snes.c
[0]PETSC ERROR: #5 main() line 808 in
/home/fuentes/github/DMPlexApplications/ireSolver.c
[0]PETSC ERROR: ----------------End of Error Message -------send entire
error message to petsc-maint at mcs.anl.gov----------
application called MPI_Abort(MPI_COMM_WORLD, 72) - process 0




thanks,
David
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140421/e7550c5c/attachment.html>


More information about the petsc-users mailing list