[petsc-users] Problem setting Jacobian matrix in complementarity problem

Barry Smith bsmith at mcs.anl.gov
Mon Apr 11 15:16:00 CDT 2016


> On Apr 11, 2016, at 2:04 PM, Christopher Athaide <cnathaide at yahoo.com> wrote:
> 
> Hi,
>  
> I am trying to solve a mixed complementarity problem using the TAO SSILS Solver . 
>  
> C(x) = A x + d ,    l <= x <= u
>  
> where A = n x n matrix
>                D = n x 1 vector
>  
> The Jacobian matrix is equal to the matrix A and is constant at each iteration.
>  
> I have setup the Jacobian update function
>  
> PetscErrorCode FormJacobian(Tao, Vec X, Mat H, Mat, void *ptr)
> {
> AppCtx         *data = (AppCtx *)ptr;
>                 PetscErrorCode ierr;
>                 PetscInt       N;
>                 PetscInt       zero = 0;
>                 PetscScalar    *d, v;
>                 PetscBool      assembled;
>  
>                 PetscFunctionBegin;
>                 ierr = MatAssembled(H, &assembled); CHKERRQ(ierr);
>                 if (assembled) { ierr = MatZeroEntries(H); CHKERRQ(ierr); }
>  
>                 /* Get pointers to vector data */
>                 VecGetSize(X, &N);
>                 for (auto i = zero; i != N; ++i)
>                 {
>                                 for (auto j = zero; j != N; ++j)
>                                 {
>                                                 ierr = MatGetValue(data->A, i, j, &v);
>                                                 ierr = MatSetValue(H, i, j, v, INSERT_VALUES);
>                                 }
>                 }
>  
>                 /* Assemble the matrix */
>                 ierr = MatSetOption(H, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE); CHKERRQ(ierr);
>                 ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
>                 ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
>                 PetscFunctionReturn(0);
>  
> }
>  
> In this case, after iterating over all the row and column indexes and assembling the matrix the function seems to run correctly.
>  
>  
> But when I try to copy the matrix
>  
> PetscErrorCode FormJacobian(Tao, Vec X, Mat H, Mat, void *ptr)
> {
>  
> PetscFunctionReturn(0);
> AppCtx         *data = (AppCtx *)ptr;
> PetscErrorCode ierr;
> ierr = MatAssembled(H, &assembled); CHKERRQ(ierr);
> if (assembled) { ierr = MatZeroEntries(H); CHKERRQ(ierr); }
> ierr = MatConvert(data->A, MATSAME, MatReuse::MAT_INITIAL_MATRIX, &H); CHKERRQ(ierr);

   You cannot do this because you are creating an entirely new H matrix here in the subroutine that TAO won't know about.

   You should use MatCopy(data->A, H, SAME_NONZERO_PATTERN);

   But since H is constant why have a copy? Just set the values in H initially and don't have a data->A matrix at all. Then your FormJacobian() doesn't need to do anything since H is already set with proper values and H remains the same through all iterations.



   Barry


>  
>                 
>  
>                 /* Assemble the matrix */
>                 ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 
>                 ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
>                 PetscFunctionReturn(0);
> }
>  
> I get the following error messages.
>  
> Problem size: 259 rows
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: Object is in wrong state
> [0]PETSC ERROR: Not for unassembled matrix
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.6.3, Dec, 03, 2015
> [0]PETSC ERROR: Unknown Name on a arch-mswin-c-debug named LAGRANGE by Christopher Mon Apr 11 19:57:37 2016
> [0]PETSC ERROR: Configure options --with-cc="win32fe cl" --with-fc=0 --with-mpi=0 --download-f2cblaslapack=1 --with-debugging=1
> [0]PETSC ERROR: #1 MatMult() line 2210 in C:\Apps\Projects\petsc-3.6.3\src\mat\interface\matrix.c
> [0]PETSC ERROR: #2 MatDFischer() line 284 in C:\Apps\Projects  \petsc-3.6.3\src\tao\util\tao_util.c
> [0]PETSC ERROR: #3 Tao_SSLS_FunctionGradient() line 64 in C:\Apps\Projects  \petsc-3.6.3\src\tao\complementarity\impls\ssls\ssls.c
> [0]PETSC ERROR: #4 TaoLineSearchComputeObjectiveAndGradient() line 941 in C:\Apps\Projects \petsc-3.6.3\src\tao\linesearch\interface\taolinesearch.c
> [0]PETSC ERROR: #5 TaoSolve_SSILS() line 64 in C:\Apps\Projects \petsc-3.6.3\src\tao\complementarity\impls\ssls\ssils.c
> [0]PETSC ERROR: #6 TaoSolve() line 204 in C:\Apps\Projects\ petsc-3.6.3\src\tao\interface\taosolver.c
> [0]PETSC ERROR: #7 run_complementarity() line 394 in c:\apps\projects \complementarity.cpp
>  
> I am new to PETSC. Can anyone help me?
>  
> Thanks
>  
> Chris Athiade



More information about the petsc-users mailing list