[petsc-users] Problem with GAMG elasticity case:'createdefaultdata' not set(?) need to support NULL data!

Matthew Knepley knepley at gmail.com
Wed Oct 16 15:07:45 CDT 2013


On Wed, Oct 16, 2013 at 3:03 PM, Mark F. Adams <mfadams at lbl.gov> wrote:

> Damn Matt, you know my code better than I do.
>
> Do we require that users call SetFromOptions? Or PCGAMGSetType in this
> case. This seems a little bad.
>

I would normally say that it should have a default GAMGType on creation.

    Matt


> Mark
>
> On Oct 16, 2013, at 2:28 PM, Matthew Knepley <knepley at gmail.com> wrote:
>
> On Wed, Oct 16, 2013 at 1:22 PM, Einar Sørheim <einar.sorheim at gmail.com>wrote:
>
>> Hi
>> I am testing out the petsc ksp solvers in our own FEM code, in particluar
>> we are interseted in the GAMG preconditioner.
>> Applying it to our mechanical problem we get some error messages from
>> petsc, which I suspect is due to some necessary init calls that are
>> missing. First the error message:
>> [0]PETSC ERROR: --------------------- Error Message
>> ------------------------------------
>> [0]PETSC ERROR: Petsc has generated inconsistent data!
>> [0]PETSC ERROR: 'createdefaultdata' not set(?) need to support NULL data!
>> [0]PETSC ERROR:
>> ------------------------------------------------------------------------
>> [0]PETSC ERROR: Petsc Release Version 3.4.3, Oct, 15, 2013
>> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
>> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
>> [0]PETSC ERROR: See docs/index.html for manual pages.
>> [0]PETSC ERROR:
>> ------------------------------------------------------------------------
>> [0]PETSC ERROR: D:\Programming\gits\Alsim\SharedModules -
>> Copy\Source\Release\Alsim_nypetsc.exe on a arch-mswin-c-opt named CMP3 by
>> admeinar Wed Oct 16 18:23:32 2013
>> [0]PETSC ERROR: Libraries linked from
>> /cygdrive/d/Programming/petsc-3.4.3/arch-mswin-c-opt/lib
>> [0]PETSC ERROR: Configure run at Wed Oct 16 14:46:01 2013
>> [0]PETSC ERROR: Configure options --with-debugging=0 --with-openmp=yes
>> --with-pthread=no --with-cc="win32fe icl" --with-fc="win32fe ifort"
>> --with-cxx="win32fe icl" --with-blas-lapack-dir="C:/Program Files
>> (x86)/Intel/Composer XE 2013/mkl/lib/intel64" --download-superlu
>> --with-sowing=0 --with-c2html=0 --with-mpi-dir="C:/Program Files
>> (x86)/Intel/MPI/4.1.0.028/em64t" --useThreads=0 --useThreads=0
>> [0]PETSC ERROR:
>> ------------------------------------------------------------------------
>> [0]PETSC ERROR: PCSetUp_GAMG() line 609 in
>> src/ksp/pc/impls/gamg/D:\PROGRA~1\PETSC-~1.3\src\ksp\pc\impls\gamg\gamg.c
>> [0]PETSC ERROR: PCSetUp() line 890 in
>> src/ksp/pc/interface/D:\PROGRA~1\PETSC-~1.3\src\ksp\pc\INTERF~1\precon.c
>> [0]PETSC ERROR: KSPSetUp() line 278 in
>> src/ksp/ksp/interface/D:\PROGRA~1\PETSC-~1.3\src\ksp\ksp\INTERF~1\itfunc.c
>>
>> In short our code looks like the following:
>>
>>
>>       call KSPGetPC(ksp,pc,ierr)
>>       tol = 1.e-20
>>       call
>> KSPSetTolerances(ksp,tol,PETSC_DEFAULT_DOUBLE_PRECISION,PETSC_DEFAULT_DOUBLE_PRECISION,imaxiter,ierr)
>>       select case (i_Precond)
>>
>>       case ( 3 )
>>          call PCSetType(pc,PCGAMG,ierr)
>> !         call PCGAMGInitializePackage()
>> !         call PCCreate_GAMG(pc,ierr)
>>       case DEFAULT
>>          call PCSetType(pc,PCJACOBI,ierr)
>>       end select
>>
>> !      call PCSetType(pc,PCJACOBI,ierr)
>>       select case (i_solver)
>>       case ( 0 )
>>          call KSPSetType(ksp,KSPCG,ierr)
>>       case DEFAULT
>>          call KSPSetType(ksp,KSPCG,ierr)
>>       end select
>>       call KSPCGSetType(ksp,KSP_CG_SYMMETRIC,ierr)
>>       if (i_nodedof==3) then
>>          write(*,*) "Set coord, number of nodes is:", i_nodes/i_nodedof
>>          call
>> VecCreateSeqWithArray(MPI_COMM_SELF,3,i_nodes,coords_,vec_coords,ierr)
>>          call MatNullSpaceCreateRigidBody(vec_coords,Matnull,ierr)
>>          call MatSetNearNullSpace(A,Matnull,ierr)
>>          call MatNullSpaceDestroy(Matnull,ierr)
>>          call PetscViewerASCIIOpen(PETSC_COMM_WORLD, "Coordvec.m",
>> viewer,ierr)
>>          call PetscViewerSetFormat( viewer,
>> PETSC_VIEWER_ASCII_MATLAB,ierr)
>>          call VecView(vec_coords,viewer,ierr)
>>          call PetscViewerDestroy( viewer,ierr )
>>          call VecDestroy(vec_coords,ierr)
>>          write(*,*) "Finish setting null space ierr =", ierr
>> !         call PCSetCoordinates( pc, 3,i_nodes/i_nodedof, coords_, ierr )
>> !         CHKERRQ(ierr)
>>       end if
>>
>>
>> !      call KSPSetFromOptions(ksp,ierr)
>>
>
> It looks like you have not set the type of aggregation. If you uncomment
> this SetFromOptions()
> it should be done automatically.
>
>    Matt
>
>
>>       call KSPSetInitialGuessNonzero(ksp,PETSC_TRUE,ierr)
>>       call KSPSetUp(ksp,ierr)
>>       call PetscTime(t0,ierr)
>>       write(*,*) "Time setup", t0-t1
>> ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
>> !                      Solve the linear system
>> ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
>>
>>       call KSPSolve(ksp,b,x,ierr)
>>       call PetscTime(t1,ierr)
>>
>> --
>> Einar Sørheim
>>
>
>
>
> --
> 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/20131016/7452d4ac/attachment.html>


More information about the petsc-users mailing list