[petsc-users] Nonconforming object sizes using TAO (TAOBNTR)

Matthew Knepley knepley at gmail.com
Fri Jan 13 19:38:53 CST 2023


On Fri, Jan 13, 2023 at 3:21 PM Alexis Marboeuf <alexis.marboeuf at hotmail.fr>
wrote:

> Hi Matt,
>
> Here is the output from Petsc when I view the TAO solver:
>

Okay, there is no sophisticated caching going on. So, first I would get it
to fail on1 process. It should if it just depends on
the convergence (I hope). Then send the source so we can run it. It should
be simple for us to find where system size
changes and the KSP is not reset (if that indeed is what happens).

  Thanks,

    Matt


> Tao Object: (Damage_) 4 MPI processes
>
>   type: bntr
>
>     Tao Object: (Damage_tao_bnk_cg_) 4 MPI processes
>
>       type: bncg
>
>         CG Type: ssml_bfgs
>
>         Skipped Stepdirection Updates: 0
>
>         Scaled gradient steps: 0
>
>         Pure gradient steps: 0
>
>         Not a descent direction: 0
>
>         Line search fails: 0
>
>         Matrix has not been preallocated yet
>
>       TaoLineSearch Object: (Damage_tao_bnk_cg_) 4 MPI processes
>
>         type: more-thuente
>
>         maximum function evaluations=30
>
>         tolerances: ftol=0.0001, rtol=1e-10, gtol=0.9
>
>         total number of function evaluations=0
>
>         total number of gradient evaluations=0
>
>         total number of function/gradient evaluations=0
>
>         Termination reason: 0
>
>       Active Set subset type: subvec
>
>       convergence tolerances: gatol=1e-08,       steptol=0.,
> gttol=0.
>
>       Residual in Function/Gradient:=0.
>
>       Objective value=0.
>
>       total number of iterations=0,                              (max:
> 2000)
>
>       Solver terminated: 0
>
>     Rejected BFGS updates: 0
>
>     CG steps: 0
>
>     Newton steps: 11
>
>     BFGS steps: 0
>
>     Scaled gradient steps: 0
>
>     Gradient steps: 0
>
>     KSP termination reasons:
>
>       atol: 4
>
>       rtol: 0
>
>       ctol: 7
>
>       negc: 0
>
>       dtol: 0
>
>       iter: 0
>
>       othr: 0
>
>   TaoLineSearch Object: (Damage_) 4 MPI processes
>
>     type: more-thuente
>
>     maximum function evaluations=30
>
>     tolerances: ftol=0.0001, rtol=1e-10, gtol=0.9
>
>     total number of function evaluations=0
>
>     total number of gradient evaluations=0
>
>     total number of function/gradient evaluations=0
>
>     using variable bounds
>
>     Termination reason: 0
>
>   KSP Object: (Damage_tao_bnk_) 4 MPI processes
>
>     type: stcg
>
>     maximum iterations=10000, nonzero initial guess
>
>     tolerances:  relative=1e-08, absolute=1e-08, divergence=1e+10
>
>     left preconditioning
>
>     using UNPRECONDITIONED norm type for convergence test
>
>   PC Object: (Damage_tao_bnk_) 4 MPI processes
>
>     type: lmvm
>
>     Mat Object: (Damage_tao_bnk_pc_lmvm_) 4 MPI processes
>
>       type: lmvmbfgs
>
>       rows=30634, cols=30634
>
>         Scale type: DIAGONAL
>
>         Scale history: 1
>
>         Scale params: alpha=1., beta=0.5, rho=1.
>
>         Convex factors: phi=0., theta=0.125
>
>         Max. storage: 5
>
>         Used storage: 5
>
>         Number of updates: 11
>
>         Number of rejects: 0
>
>         Number of resets: 0
>
>         Mat Object: (Damage_tao_bnk_pc_lmvm_J0_) 4 MPI processes
>
>           type: lmvmdiagbroyden
>
>           rows=30634, cols=30634
>
>             Scale history: 1
>
>             Scale params: alpha=1., beta=0.5, rho=1.
>
>             Convex factor: theta=0.125
>
>             Max. storage: 1
>
>             Used storage: 1
>
>             Number of updates: 11
>
>             Number of rejects: 0
>
>             Number of resets: 0
>
>     linear system matrix = precond matrix:
>
>     Mat Object: 4 MPI processes
>
>       type: mpiaij
>
>       rows=468, cols=468
>
>       total: nonzeros=2932, allocated nonzeros=2932
>
>       total number of mallocs used during MatSetValues calls=0
>
>         not using I-node (on process 0) routines
>
>   total KSP iterations: 103
>
>   Active Set subset type: subvec
>
>   convergence tolerances: gatol=0.0001,   steptol=0.,   gttol=1e-05
>
>   Residual in Function/Gradient:=9.11153e-05
>
>   Objective value=0.00665458
>
>   total number of iterations=11,                          (max: 50)
>
>   total number of function evaluations=17,                  max: -1
>
>   total number of gradient evaluations=13,                  max: -1
>
>   total number of Hessian evaluations=12
>
>   Solution converged:    ||g(X)|| <= gatol
>
> Thanks again for your help!
> Alexis
>
> ------------------------------
> *De :* Matthew Knepley <knepley at gmail.com>
> *Envoyé :* samedi 14 janvier 2023 01:38
> *À :* Alexis Marboeuf <alexis.marboeuf at hotmail.fr>
> *Cc :* petsc-users at mcs.anl.gov <petsc-users at mcs.anl.gov>
> *Objet :* Re: [petsc-users] Nonconforming object sizes using TAO (TAOBNTR)
>
> On Fri, Jan 13, 2023 at 11:22 AM Alexis Marboeuf <
> alexis.marboeuf at hotmail.fr> wrote:
>
> Hi all,
>
> In a variational approach of brittle fracture setting, I try to solve a
> bound constraint minimization problem using TAO. I checkout on the main
> branch of Petsc. Minimization with respect to the bounded variable (damage)
> is achieved through the Bounded Newton Trust Region (TAOBNTR). All other
> TAO parameters are set by default. On a Linux machine, I get the following
> error with a 4 processors run:
>
>
> Can you view the solver?
>
>   Thanks,
>
>     Matt
>
>
>
> [3]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
> [3]PETSC ERROR: Nonconforming object sizes
> [3]PETSC ERROR: Preconditioner number of local rows 1122 does not equal
> input vector size 1161
> [3]PETSC ERROR: See https://petsc.org/release/faq/ for trouble shooting.
> [3]PETSC ERROR: Petsc Development GIT revision: v3.18.3-342-gdab44c92d91
> GIT Date: 2023-01-04 13:37:04 +0000
> [3]PETSC ERROR:
> /home/marboeua/Developpement/mef90/arch-darwin-c/bin/vDefTAO on a
> arch-darwin-c named bb01 by marboeua Thu Jan 12 16:55:18 2023
> [2]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
> [3]PETSC ERROR: Configure options --FFLAGS=-ffree-line-length-none
> --COPTFLAGS="-O3 -march=znver3 -g" --CXXOPTFLAGS="-O3 -march=znver3 -g"
>  --FOPTFLAGS="-O3 -march=znver3 -g" --download-fblaslapack=1
>  --download-mumps=1 --download-chaco=1 --download-exodusii=1
>  --download-hypre=1 --download-ml=1 --download-triangle
> --download-scalapack=1 --download-superlu=1 --download-sowing=1
>  --download-sowing-cc=/opt/rh/devtoolset-9/root/usr/bin/gcc
> --download-sowing-cxx=/opt/rh/devtoolset-9/root/usr/bin/g++
> --download-sowing-cpp=/opt/rh/devtoolset-9/root/usr/bin/cpp
> --download-sowing-cxxcpp=/opt/rh/devtoolset-9/root/usr/bin/cpp
> --download-yaml=1 --download-bison=1 --download-hdf5=1 --download-metis=1
>  --download-parmetis=1 --download-netcdf=1 --download-pnetcdf=1
>  --download-zlib=1 --with-cmake=1 --with-debugging=0
>  --with-mpi-dir=/opt/HPC/mvapich2/2.3.7-gcc11.2.1 --with-ranlib=ranlib
> --with-shared-libraries=1 --with-sieve=1 --download-p4est=1 --with-pic
> --with-mpiexec=srun --with-x11=0 PETSC_ARCH=arch-darwin-c
> [3]PETSC ERROR: #1 PCApply() at /1
> /home/marboeua/Developpement/petsc/src/ksp/pc/interface/precon.c:434
> [3]PETSC ERROR: #2 KSP_PCApply() at
> /home/marboeua/Developpement/petsc/include/petsc/private/kspimpl.h:380
> [3]PETSC ERROR: #3 KSPCGSolve_STCG() at /1
> /home/marboeua/Developpement/petsc/src/ksp/ksp/impls/cg/stcg/stcg.c:76
> [3]PETSC ERROR: #4 KSPSolve_Private() at /1
> /home/marboeua/Developpement/petsc/src/ksp/ksp/interface/itfunc.c:898
> [3]PETSC ERROR: #5 KSPSolve() at /1
> /home/marboeua/Developpement/petsc/src/ksp/ksp/interface/itfunc.c:1070
> [3]PETSC ERROR: #6 TaoBNKComputeStep() at /1
> /home/marboeua/Developpement/petsc/src/tao/bound/impls/bnk/bnk.c:459
> [3]PETSC ERROR: #7 TaoSolve_BNTR() at /1
> /home/marboeua/Developpement/petsc/src/tao/bound/impls/bnk/bntr.c:138
> [3]PETSC ERROR: #8 TaoSolve() at /1
> /home/marboeua/Developpement/petsc/src/tao/interface/taosolver.c:177
> [2]PETSC ERROR: Nonconforming object sizes
> [2]PETSC ERROR: Preconditioner number of local rows 1229 does not equal
> input vector size 1254
> [2]PETSC ERROR: See https://petsc.org/release/faq/ for trouble shooting.
> [2]PETSC ERROR: Petsc Development GIT revision: v3.18.3-342-gdab44c92d91
> GIT Date: 2023-01-04 13:37:04 +0000
> [2]PETSC ERROR:
> /home/marboeua/Developpement/mef90/arch-darwin-c/bin/vDefTAO on a
> arch-darwin-c named bb01 by marboeua Thu Jan 12 16:55:18 2023
> [2]PETSC ERROR: Configure options --FFLAGS=-ffree-line-length-none
> --COPTFLAGS="-O3 -march=znver3 -g" --CXXOPTFLAGS="-O3 -march=znver3 -g"
>  --FOPTFLAGS="-O3 -march=znver3 -g" --download-fblaslapack=1
>  --download-mumps=1 --download-chaco=1 --download-exodusii=1
>  --download-hypre=1 --download-ml=1 --download-triangle
> --download-scalapack=1 --download-superlu=1 --download-sowing=1
>  --download-sowing-cc=/opt/rh/devtoolset-9/root/usr/bin/gcc
> --download-sowing-cxx=/opt/rh/devtoolset-9/root/usr/bin/g++
> --download-sowing-cpp=/opt/rh/devtoolset-9/root/usr/bin/cpp
> --download-sowing-cxxcpp=/opt/rh/devtoolset-9/root/usr/bin/cpp
> --download-yaml=1 --download-bison=1 --download-hdf5=1 --download-metis=1
>  --download-parmetis=1 --download-netcdf=1 --download-pnetcdf=1
>  --download-zlib=1 --with-cmake=1 --with-debugging=0
>  --with-mpi-dir=/opt/HPC/mvapich2/2.3.7-gcc11.2.1 --with-ranlib=ranlib
> --with-shared-libraries=1 --with-sieve=1 --download-p4est=1 --with-pic
> --with-mpiexec=srun --with-x11=0 PETSC_ARCH=arch-darwin-c
> [2]PETSC ERROR: #1 PCApply() at /1
> /home/marboeua/Developpement/petsc/src/ksp/pc/interface/precon.c:434
> [2]PETSC ERROR: #2 KSP_PCApply() at
> /home/marboeua/Developpement/petsc/include/petsc/private/kspimpl.h:380
> [2]PETSC ERROR: #3 KSPCGSolve_STCG() at /1
> /home/marboeua/Developpement/petsc/src/ksp/ksp/impls/cg/stcg/stcg.c:76
> [2]PETSC ERROR: #4 KSPSolve_Private() at /1
> /home/marboeua/Developpement/petsc/src/ksp/ksp/interface/itfunc.c:898
> [2]PETSC ERROR: #5 KSPSolve() at /1
> /home/marboeua/Developpement/petsc/src/ksp/ksp/interface/itfunc.c:1070
> [2]PETSC ERROR: #6 TaoBNKComputeStep() at /1
> /home/marboeua/Developpement/petsc/src/tao/bound/impls/bnk/bnk.c:459
> [2]PETSC ERROR: #7 TaoSolve_BNTR() at /1
> /home/marboeua/Developpement/petsc/src/tao/bound/impls/bnk/bntr.c:138
> [2]PETSC ERROR: #8 TaoSolve() at /1
> /home/marboeua/Developpement/petsc/src/tao/interface/taosolver.c:177
> [3]PETSC ERROR: #9 /home/marboeua/Developpement/mef90/vDef/vDefTAO.F90:370
> application called MPI_Abort(MPI_COMM_SELF, 60) - process 0
> [2]PETSC ERROR: #9 /home/marboeua/Developpement/mef90/vDef/vDefTAO.F90:370
> application called MPI_Abort(MPI_COMM_SELF, 60) - process 0
> slurmstepd: error: *** STEP 5034.0 ON bb01 CANCELLED AT
> 2023-01-12T17:21:07 ***
> srun: Job step aborted: Waiting up to 32 seconds for job step to finish.
> srun: error: bb01: tasks 0-1: Killed
> srun: error: bb01: tasks 2-3: Exited with exit code 1
>
> The error is raised in the middle of the computation after many successful
> calls of TAOSolve and TAO iterations. My guess is that TAO computes the
> preconditioner during its first iteration with all variables in the active
> set. But the preconditioner is never updated when some variables are moved
> to the inactive set during the next TAO iterations. Am I right? Can you
> help me with that?
>
> Thanks a lot for your help and your time.
> Regards,
> Alexis
>
>
>
> --
> 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
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>


-- 
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

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230113/3034f3b5/attachment-0001.html>


More information about the petsc-users mailing list