[petsc-users] Copy stiffness matrix to jacobian matrix in form jacobian

Barry Smith bsmith at petsc.dev
Mon Jun 14 20:04:25 CDT 2021


  How are you obtaining the original form of the copied matrix? With MatDuplicate(), do you use MAT_DO_NOT_COPY_VALUES? When you do the MatCopy() what do you pass for the MatStructure argument? SAME_NONZERO_PATTERN?  It seems like the options you are using is causing the zero locations to not be copied over.

   You are using a rather old PETSc, it should work ok but it is always best, when possible to develop code using the newest release.

  Barry


> On Jun 14, 2021, at 4:13 PM, Kaushik Vijaykumar <kaushikv318 at gmail.com> wrote:
> 
> Thanks Matt and Barry for your response. I did not copy paste the entire error as it was very long, please see the error below.
> 
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: Object is in wrong state
> [0]PETSC ERROR: Matrix is missing diagonal entry 0
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html <http://www.mcs.anl.gov/petsc/documentation/faq.html> for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.11.2, May, 18, 2019 
> [0]PETSC ERROR: ./openspaf090.x on a  named nid00060 by de764124 Mon Jun 14 07:18:56 2021
> [0]PETSC ERROR: Configure options --known-has-attribute-aligned=1 --known-mpi-int64_t=0 --known-bits-per-byte=8 --known-64-bit-blas-indices=0 --known-sdot-returns-double=0 --known-snrm2-returns-double=0 --known-level1-dcache-assoc=0 --known-level1-dcache-linesize=32 --known-level1-dcache-size=32768 --known-memcmp-ok=1 --known-mpi-c-double-complex=1 --known-mpi-long-double=1 --known-mpi-shared-libraries=0 --known-sizeof-MPI_Comm=4 --known-sizeof-MPI_Fint=4 --known-sizeof-char=1 --known-sizeof-double=8 --known-sizeof-float=4 --known-sizeof-int=4 --known-sizeof-long-long=8 --known-sizeof-long=8 --known-sizeof-short=2 --known-sizeof-size_t=8 --known-sizeof-void-p=8 --with-ar=ar --with-batch=1 --with-cc=cc --with-clib-autodetect=0 --with-cxx=CC --with-cxx-dialect=C++11 --with-cxxlib-autodetect=0 --with-debugging=0 --with-dependencies=0 --with-fc=ftn --with-fortran-datatypes=0 --with-fortran-interfaces=0 --with-fortranlib-autodetect=0 --with-ranlib=ranlib --with-scalar-type=real --with-shared-ld=ar --with-etags=0 --with-dependencies=0 --with-x=0 --with-ssl=0 --with-shared-libraries=0 --with-dependencies=0 --with-mpi-lib="[]" --with-mpi-include="[]" --with-blaslapack-lib="-L/opt/cray/pe/libsci/19.04.1.1/INTEL/16.0/x86_64/lib <http://19.4.1.1/INTEL/16.0/x86_64/lib> -lsci_intel_mp" --with-superlu=1 --with-superlu-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-superlu-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lsuperlu" --with-superlu_dist=1 --with-superlu_dist-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-superlu_dist-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lsuperlu_dist" --with-parmetis=1 --with-parmetis-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-parmetis-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lparmetis" --with-metis=1 --with-metis-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-metis-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lmetis" --with-ptscotch=1 --with-ptscotch-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-ptscotch-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lptscotch -lscotch -lptscotcherr -lscotcherr" --with-scalapack=1 --with-scalapack-include=/opt/cray/pe/libsci/19.04.1.1/INTEL/16.0/x86_64/include <http://19.4.1.1/INTEL/16.0/x86_64/include> --with-scalapack-lib="-L/opt/cray/pe/libsci/19.04.1.1/INTEL/16.0/x86_64/lib <http://19.4.1.1/INTEL/16.0/x86_64/lib> -lsci_intel_mpi_mp -lsci_intel_mp" --with-mumps=1 --with-mumps-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-mumps-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -L/opt/cray/pe/mpt/7.7.7/gni/mpich-intel/16.0/lib -lcmumps -ldmumps -lesmumps -lsmumps -lzmumps -lmumps_common -lptesmumps -lesmumps -lpord -lfmpich" --with-hdf5=1 --with-hdf5-include=/opt/cray/pe/hdf5-parallel/1.10.5.0/intel/18.0/include <http://1.10.5.0/intel/18.0/include> --with-hdf5-lib="-L/opt/cray/pe/hdf5-parallel/1.10.5.0/intel/18.0/lib <http://1.10.5.0/intel/18.0/lib> -lhdf5_parallel -lz -ldl" --CFLAGS="-xavx -qopenmp -O3  -fpic" --CPPFLAGS= --CXXFLAGS="-xavx -qopenmp -O3   -fpic" --FFLAGS="-xavx -qopenmp -O3  -fpic" --LIBS=-lstdc++ --CXX_LINKER_FLAGS= --PETSC_ARCH=sandybridge --prefix=/opt/cray/pe/petsc/3.11.2.0/real/INTEL/19.0/sandybridge <http://3.11.2.0/real/INTEL/19.0/sandybridge> --with-hypre=1 --with-hypre-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-hypre-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lHYPRE" --with-sundials=1 --with-sundials-include=/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/include --with-sundials-lib="-L/opt/cray/pe/tpsl/19.06.1/INTEL/19.0/sandybridge/lib -lsundials_cvode -lsundials_cvodes -lsundials_ida -lsundials_idas -lsundials_kinsol -lsundials_nvecparallel -lsundials_nvecserial"
> [0]PETSC ERROR: #1 MatILUFactorSymbolic_SeqAIJ() line 1687 in src/mat/impls/aij/seq/aijfact.c
> [0]PETSC ERROR: #2 MatILUFactorSymbolic() line 6759 in src/mat/interface/matrix.c
> [0]PETSC ERROR: #3 PCSetUp_ILU() line 144 in src/ksp/pc/impls/factor/ilu/ilu.c
> [0]PETSC ERROR: #4 PCSetUp() line 932 in src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #5 KSPSetUp() line 391 in src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #6 KSPSolve() line 725 in src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #7 SNESSolve_NEWTONLS() line 225 in src/snes/impls/ls/ls.c
> [0]PETSC ERROR: #8 SNESSolve() line 4560 in src/snes/interface/snes.c
> [0]PETSC ERROR: #9 main() line 601 in spaf090_main.c
> [0]PETSC ERROR: No PETSc Option Table entries
> [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
> Rank 0 [Mon Jun 14 07:18:57 2021] [c0-0c0s15n0] application called MPI_Abort(MPI_COMM_WORLD, 73) - process 0
> _pmiu_daemon(SIGCHLD): [NID 00060] [c0-0c0s15n0] [Mon Jun 14 07:18:57 2021] PE RANK 0 exit signal Aborted
> Application 11927 exit codes: 134
> Application 11927 resources: utime ~0s, stime ~1s, Rss ~783564, inblocks ~0, outblocks ~0
> 
> 
> In my previous approach, I was copying the jacobian in formjacobian and I did not encounter this error, of Matrix is missing a diagonal entry 0
> 
> If I copy my jacobian using MatGetValue and MatSetValue, the above error disappears
> 
> formjacobian()
> {
> 
> my_ctx* ptr = (my_ctx*) ctx; // cast the pointer to void into pointer to struct
> 
> ierr = MatGetOwnershipRange(jac,&istart,&iend);
> ierr = MatGetSize(ptr->K,&m,&n); CHKERRQ(ierr);
> 
> ierr = MatAssemblyBegin(ptr->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> ierr = MatAssemblyEnd(ptr->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> 
>  for(i=istart;i<iend;i++)
> {
>  for(j=0;j<n;j++)
>  {
>  ierr = MatGetValue(ptr->K, i, j, &v);CHKERRQ(ierr);
>  ierr = MatSetValue(jac, i, j, v, INSERT_VALUES);CHKERRQ(ierr);
>  ierr = MatSetValue(B, i, j, v, INSERT_VALUES);CHKERRQ(ierr);
>  }
>  }
> 
> ...
> }
> 
> 
> Thanks
> Kaushik
> 
> On Mon, Jun 14, 2021 at 3:45 PM Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>> wrote:
> 
>   It is best to cut and paste the entire error message, not just two lines. We are missing all context to see where the error is occurring and possible causes.
> 
> 
> 
>> On Jun 14, 2021, at 6:20 AM, Kaushik Vijaykumar <kaushikv318 at gmail.com <mailto:kaushikv318 at gmail.com>> wrote:
>> 
>> Thanks for the suggestion Barry, I did do that by setting it in multiple ways first by setting
>> 
>> Formfunction(snes,Vec x,Vec F,void* ctx)
>> {
>> ierr = SNESGetJacobian(snes,&jac,NULL,NULL,NULL);CHKERRQ(ierr);
>> ierr = MatZeroEntries(jac);CHKERRQ(ierr);
>> 
>> .....
>> 
>> 
>> ierr = MatGetOwnershipRange(jac,&istart,&iend);
>> 
>> for(i=istart;i<iend;i++)
>> {
>> val = 0.0;
>> ierr = MatSetValue(jac, i, i, val, ADD_VALUES);CHKERRQ(ierr);
>> }
>> ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>> ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>> }
>> I get the same error as before for a single processor run.
>> "[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
>> [0]PETSC ERROR: Object is in wrong state
>> [0]PETSC ERROR: Matrix is missing diagonal entry 0"
>> 
>> 
>> Thanks
>> Kaushik
>> 
>> On Fri, Jun 11, 2021 at 11:33 AM Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>> wrote:
>> 
>>> When I execute this I get a Petsc error, object in wrong state missing diagonal entry.
>> 
>>    When do you get this error? During a TSSolve? While KSP is building a preconditioner.
>> 
>>    Some parts of PETSc require that you put entries (even if they are zero) on all diagonal entries in the Jacobian. So likely you need to just make sure that your populate jac routine puts/adds 0 to all diagonal entries.
>> 
>>   Barry
>> 
>> 
>> 
>>> On Jun 11, 2021, at 10:16 AM, Kaushik Vijaykumar <kaushikv318 at gmail.com <mailto:kaushikv318 at gmail.com>> wrote:
>>> 
>>> Thanks Barry for your comments. Based on your suggestion I tried the following
>>> 
>>> Main()
>>> {
>>> ierr = SNESSetFunction(snes,r1,FormFunction,&this_ctx);CHKERRQ(ierr);
>>> ierr = SNESSetJacobian(snes,jac,NULL,FormJacobian,&this_ctx);CHKERRQ(ierr);
>>> ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
>>> ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
>>> }
>>> 
>>> In formfunction I populate the jacobian "jac"
>>> 
>>> Formfunction(snes,Vec x,Vec F,void* ctx)
>>> {
>>> ierr = SNESGetJacobian(snes,&jac,Null,NULL,NULL);CHKERRQ(ierr);
>>> 
>>> // "Populate jac"
>>> ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>>> ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);  
>>> 
>>> }
>>> FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
>>> {
>>> ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>>> ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>>> 
>>> ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>>> ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>>> }
>>> 
>>> When I execute this I get a Petsc error, object in wrong state missing diagonal entry.
>>> 
>>> I am not sure what I am missing here? 
>>> 
>>> Thanks
>>> Kaushik
>>> 
>>> 
>>> On Thu, Jun 10, 2021 at 6:37 PM Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>> wrote:
>>> 
>>> 
>>> > On Jun 10, 2021, at 3:29 PM, Kaushik Vijaykumar <kaushikv318 at gmail.com <mailto:kaushikv318 at gmail.com>> wrote:
>>> > 
>>> > Hi everyone,
>>> > 
>>> > I am trying to copy the stiffness matrix that I generated in form function to the jacobian matrix jac in form jacobian. The piece of code for that:
>>> > 
>>> > 
>>>    If you have created jac and B in your main program and passed them with SNESSetJacobian() then you should be able to just use
>>> 
>>> > ierr = MatCopy(ptr->K,jac,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
>>> 
>>> Generally jac and B are the same matrix so you don't need the second copy. 
>>> > ierr = MatCopy(ptr->K,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);  
>>> 
>>> If they are sometimes the same and sometime not you can do
>>> 
>>> if (jac != B) {ierr = MatCopy(ptr->K,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);  }
>>> 
>>> 
>>>    The following won't work because you are creating new local matrices in your form jacobian that don't exist outside of that routine.
>>> 
>>> > 
>>> > ierr = MatDuplicate(ptr->K,MAT_COPY_VALUES,&jac);CHKERRQ(ierr);
>>> > ierr = MatDuplicate(ptr->K,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
>>> > 
>>> > ierr = MatCopy(ptr->K,jac,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
>>> > ierr = MatCopy(ptr->K,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);  
>>> > 
>>> > 
>>>    If the matrix you are creating in form function is the Jacobian then you don't need to use the memory and time to make copies. You can just form the matrix directly into the Jacobian you will use for SNES. You can obtain this with a call 
>>> 
>>>    SNESGetJacobian(snes,&jac,NULL,NULL,NULL); 
>>> 
>>> in your form function and just put the matrix values directly into jac. Then your form Jacobian does not need to do anything but return since the jac already contains the Jacobian.
>>> 
>>> On the off-change you are using the iteration A(x^n) x^{n+1} = b(x^n} and not using Newtons' method (which is why you compute A in form function, then you might be better off using SNESSetPicard() instead of SNESSetFunction(), SNESSetJacobian().
>>> 
>>> > 
>>> > Thanks
>>> > Kaushik
>>> 
>> 
> 

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20210614/5b492843/attachment-0001.html>


More information about the petsc-users mailing list