[petsc-users] Problem with TS and PCMG

Torquil Macdonald Sørensen torquil at gmail.com
Wed Dec 18 15:11:53 CST 2013


On 18/12/13 21:19, Jed Brown wrote:
>> In short, I'm wondering if it is possible to use TSSolve without using
>> any DM structure:
>>
>> When I was debugging my program I tried to analyse the difference between
>> the TSSolve-based time iteration program, and a version where I do the
>> time-stepping manually using a for-loop which repeateed KSPSolve functions.
>> The debugger showed that in the TSSolve program a pointer named
>> ts->snes->ksp->pc->dm was nonzero, but in the for-loop program it was zero.
>> In the TSSolve-program, and pc->dm != 0 made it enter a block of code
>> within PCSetUp_MG that called DMCoarsen, which then crashed. So I guess I'm
>> wondering how to tell PETSc that I'm not using any DM structure.
> Please show the entire trace.  Callback pointers are stored on an
> "internal" DM if you don't provide one.  This is important for nonlinear
> preconditioning and other types of composition.  What preconditioner are
> you intending to use?

If you mean the entire backtrace within GDB, I'm including it in this
email. I'm using PCMG. My minimal testcase program has everything
hardcoded. I'm not providing any runtime options other than
-start_in_debugger, so everything can be seen from the code. I've since
made a couple of small changes that might influence line number
statements from GDB, so I'm including the newest version in this message.

Here is the GDB output of "bt full" and the corresponding program code:

(gdb) bt full
#1  0x00007ffdc080d172 in DMCoarsen (dm=0x1fa7df0, comm=0x7ffdbffdd9c0
<ompi_mpi_comm_null>, dmc=0x2048f00) at dm.c:1811
        ierr = 0
        link = 0x7ffdc1af1640
        __func__ = "DMCoarsen"
#2  0x00007ffdc0c41c1e in PCSetUp_MG (pc=0x1ff7570) at mg.c:593
        kdm = 0x1d5b370
        p = 0xb453
        rscale = 0x1ff7570
        lu = PETSC_FALSE
        svd = (PETSC_TRUE | unknown: 32764)
        dB = 0x0
        __func__ = "PCSetUp_MG"
        tvec = 0x0
        viewer = 0x0
        mglevels = 0x2006980
        ierr = 0
        preonly = PETSC_FALSE
        cholesky = (unknown: 3246112504)
        mg = 0x20060c0
        i = 0
        n = 2
        redundant = (unknown: 1211234)
        use_amat = PETSC_TRUE
        dms = 0x2048f00
        cpc = 0x200f620
        dump = PETSC_FALSE
        opsset = PETSC_FALSE
        dA = 0x7ffdc1af1640
        uflag = DIFFERENT_NONZERO_PATTERN
#3  0x00007ffdc0f3b6ce in PCSetUp (pc=0x1ff7570) at precon.c:890
        ierr = 0
        def = 0x7ffdc1494e01 "User provided function"
        __func__ = "PCSetUp"
#4  0x00007ffdc0abf7e8 in KSPSetUp (ksp=0x1fce550) at itfunc.c:278
        ierr = 0
        A = 0x55
        __func__ = "KSPSetUp"
        B = 0x1fce550
        stflg = SAME_NONZERO_PATTERN
#5  0x00007ffdc0ac0a30 in KSPSolve (ksp=0x1fce550, b=0x2033dd0,
x=0x204b120) at itfunc.c:399
        rank = 0
        flag1 = PETSC_FALSE
        guess_zero = (unknown: 33861920)
        viewer = 0x7fffbf50c290
        ierr = 0
        flag2 = PETSC_FALSE
        __func__ = "KSPSolve"
        flag3 = PETSC_TRUE
        flg = PETSC_FALSE
        inXisinB = PETSC_FALSE
        mat = 0x0
        premat = 0x0
        format = PETSC_VIEWER_DEFAULT
#6  0x00007ffdc0adae0c in SNESSolve_KSPONLY (snes=0x1f86390) at ksponly.c:44
        ierr = 0
        lits = 32767
        flg = SAME_PRECONDITIONER
        Y = 0x204b120
        __func__ = "SNESSolve_KSPONLY"
        X = 0x2036b50
        F = 0x2033dd0
        kspreason = KSP_CONVERGED_ITERATING
#7  0x00007ffdc11cb642 in SNESSolve (snes=0x1f86390, b=0x0, x=0x2036b50)
at snes.c:3636
        ierr = 0
        flg = PETSC_FALSE
        viewer = 0x7ffdc18e5205 <_dl_runtime_resolve+53>
        grid = 0
        xcreated = 0x0
        __func__ = "SNESSolve"
        dm = 0x1fa7df0
        format = PETSC_VIEWER_DEFAULT
#8  0x00007ffdc12780cc in TSStep_Theta (ts=0x1f52ea0) at theta.c:183
        shift = 10
        th = 0x1fa3130
        next_time_step = 0.10000000000000001
        snesreason = 32765
        adapt = 0x7ffdc18ded95 <_dl_fixup+245>
        accept = (PETSC_TRUE | unknown: 32766)
        its = 32765
        reject = 0
        next_scheme = 0
        ierr = 0
        __func__ = "TSStep_Theta"
        lits = -1071406416
#9  0x00007ffdc12ba202 in TSStep (ts=0x1f52ea0) at ts.c:2507
        ptime_prev = 0
        ierr = 0
        __func__ = "TSStep"
#10 0x00000000004026fd in main (argc=3, argv=0x7fffbf50c6d8) at
/home/tmac/programming/petsc/source/tssolve_multigrid.cpp:65
        ts = 0x1f52ea0
        ksp = 0x1fce550
        pc = 0x1ff7570
        ierr = 0
        __func__ = "main"
        nb_dof_fine = 5
        nb_dof_coarse = 3
        x = 0x1d769d0
        Id = 0x1dffc30
        interpolation = 0x1e24020
        restriction = 0x1e48470

$ cat /home/tmac/programming/petsc/source/tssolve_multigrid.cpp
// A trivial test case for TSBEULER and PCMG: Solve Id*dU/dt = 0, giving
a segfault

#include "petscts.h"

int main(int argc, char ** argv)
{
    PetscErrorCode ierr = PetscInitialize(&argc, &argv, 0, 0);
CHKERRQ(ierr);

    int const nb_dof_fine = 5;
    int const nb_dof_coarse = 3;

    Vec x; // A vector to hold the solution on the fine grid
    ierr = VecCreateSeq(PETSC_COMM_SELF, nb_dof_fine, &x); CHKERRQ(ierr);

    Mat Id; // Identity matrix on the vector space of the fine grid
    ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nb_dof_fine, nb_dof_fine, 1,
0, &Id); CHKERRQ(ierr);
    for(int i = 0; i != nb_dof_fine; i++) { MatSetValue(Id, i, i, 1.0,
INSERT_VALUES); }
    ierr = MatAssemblyBegin(Id, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
    ierr = MatAssemblyEnd(Id, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);

    Mat interpolation; // Multigrid interpolation matrix
    ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nb_dof_fine, nb_dof_coarse,
2, 0, &interpolation); CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 0, 0, 1.0, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 1, 0, 0.5, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 1, 1, 0.5, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 2, 1, 1.0, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 3, 1, 0.5, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 3, 2, 0.5, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(interpolation, 4, 2, 1.0, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatAssemblyBegin(interpolation, MAT_FINAL_ASSEMBLY);
CHKERRQ(ierr);
    ierr = MatAssemblyEnd(interpolation, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);

    Mat restriction; // Multigrid restriction matrix, not the transpose
of "interpolation"
    ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nb_dof_coarse, nb_dof_fine,
1, 0, &restriction); CHKERRQ(ierr);
    ierr = MatSetValue(restriction, 0, 0, 1.0, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(restriction, 1, 2, 1.0, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatSetValue(restriction, 2, 4, 1.0, INSERT_VALUES);
CHKERRQ(ierr);
    ierr = MatAssemblyBegin(restriction, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
    ierr = MatAssemblyEnd(restriction, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);

    {
        TS ts;
        ierr = TSCreate(PETSC_COMM_SELF, &ts); CHKERRQ(ierr);
        ierr = TSSetProblemType(ts, TS_LINEAR); CHKERRQ(ierr);
        ierr = TSSetType(ts, TSBEULER); // Use an implicit scheme
        ierr = TSSetInitialTimeStep(ts, 0.0, 0.1); CHKERRQ(ierr);
        ierr = TSSetSolution(ts, x); CHKERRQ(ierr);
        ierr = TSSetIFunction(ts, 0, TSComputeIFunctionLinear, 0);
CHKERRQ(ierr);
        ierr = TSSetIJacobian(ts, Id, Id, TSComputeIJacobianConstant, 0);

        KSP ksp;
        ierr = TSGetKSP(ts, &ksp); CHKERRQ(ierr);
        ierr = KSPSetType(ksp, KSPFGMRES); CHKERRQ(ierr);

        PC pc;
        ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
        ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
        ierr = PCMGSetLevels(pc, 2, 0); CHKERRQ(ierr);
        ierr = PCMGSetGalerkin(pc, PETSC_TRUE); CHKERRQ(ierr);
        ierr = PCMGSetInterpolation(pc, 1, interpolation); CHKERRQ(ierr);
        ierr = PCMGSetRestriction(pc, 1, restriction); CHKERRQ(ierr);

        ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
        ierr = TSSetUp(ts); CHKERRQ(ierr);
        ierr = TSStep(ts); CHKERRQ(ierr);

        ierr = TSDestroy(&ts); CHKERRQ(ierr);
    }

    ierr = MatDestroy(&restriction); CHKERRQ(ierr);
    ierr = MatDestroy(&interpolation); CHKERRQ(ierr);
    ierr = VecDestroy(&x); CHKERRQ(ierr);
    ierr = MatDestroy(&Id); CHKERRQ(ierr);
    ierr = PetscFinalize(); CHKERRQ(ierr);

    return 0;
}

Best regards
Torquil Sørensen

-------------- next part --------------
A non-text attachment was scrubbed...
Name: signature.asc
Type: application/pgp-signature
Size: 897 bytes
Desc: OpenPGP digital signature
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20131218/13be8e7c/attachment.pgp>


More information about the petsc-users mailing list