1 2 #include <../src/snes/impls/vi/rs/virsimpl.h> /*I "petscsnes.h" I*/ 3 #include <petsc-private/kspimpl.h> 4 #include <petsc-private/matimpl.h> 5 #include <petsc-private/dmimpl.h> 6 #include <petsc-private/vecimpl.h> 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "SNESVIGetInactiveSet" 10 /* 11 SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear 12 system is solved on) 13 14 Input parameter 15 . snes - the SNES context 16 17 Output parameter 18 . ISact - active set index set 19 20 */ 21 PetscErrorCode SNESVIGetInactiveSet(SNES snes,IS *inact) 22 { 23 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*)snes->data; 24 25 PetscFunctionBegin; 26 *inact = vi->IS_inact_prev; 27 PetscFunctionReturn(0); 28 } 29 30 /* 31 Provides a wrapper to a DM to allow it to be used to generated the interpolation/restriction from the DM for the smaller matrices and vectors 32 defined by the reduced space method. 33 34 Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints. 35 36 <*/ 37 typedef struct { 38 PetscInt n; /* size of vectors in the reduced DM space */ 39 IS inactive; 40 41 PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*); /* DM's original routines */ 42 PetscErrorCode (*coarsen)(DM, MPI_Comm, DM*); 43 PetscErrorCode (*createglobalvector)(DM,Vec*); 44 45 DM dm; /* when destroying this object we need to reset the above function into the base DM */ 46 } DM_SNESVI; 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "DMCreateGlobalVector_SNESVI" 50 /* 51 DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space 52 53 */ 54 PetscErrorCode DMCreateGlobalVector_SNESVI(DM dm,Vec *vec) 55 { 56 PetscErrorCode ierr; 57 PetscContainer isnes; 58 DM_SNESVI *dmsnesvi; 59 60 PetscFunctionBegin; 61 ierr = PetscObjectQuery((PetscObject)dm,"VI",(PetscObject*)&isnes);CHKERRQ(ierr); 62 if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Composed SNES is missing"); 63 ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi);CHKERRQ(ierr); 64 ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),dmsnesvi->n,PETSC_DETERMINE,vec);CHKERRQ(ierr); 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "DMCreateInterpolation_SNESVI" 70 /* 71 DMCreateInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints. 72 73 */ 74 PetscErrorCode DMCreateInterpolation_SNESVI(DM dm1,DM dm2,Mat *mat,Vec *vec) 75 { 76 PetscErrorCode ierr; 77 PetscContainer isnes; 78 DM_SNESVI *dmsnesvi1,*dmsnesvi2; 79 Mat interp; 80 81 PetscFunctionBegin; 82 ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject*)&isnes);CHKERRQ(ierr); 83 if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm1),PETSC_ERR_PLIB,"Composed VI data structure is missing"); 84 ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr); 85 ierr = PetscObjectQuery((PetscObject)dm2,"VI",(PetscObject*)&isnes);CHKERRQ(ierr); 86 if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm2),PETSC_ERR_PLIB,"Composed VI data structure is missing"); 87 ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi2);CHKERRQ(ierr); 88 89 ierr = (*dmsnesvi1->createinterpolation)(dm1,dm2,&interp,NULL);CHKERRQ(ierr); 90 ierr = MatGetSubMatrix(interp,dmsnesvi2->inactive,dmsnesvi1->inactive,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr); 91 ierr = MatDestroy(&interp);CHKERRQ(ierr); 92 *vec = 0; 93 PetscFunctionReturn(0); 94 } 95 96 extern PetscErrorCode DMSetVI(DM,IS); 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "DMCoarsen_SNESVI" 100 /* 101 DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set 102 103 */ 104 PetscErrorCode DMCoarsen_SNESVI(DM dm1,MPI_Comm comm,DM *dm2) 105 { 106 PetscErrorCode ierr; 107 PetscContainer isnes; 108 DM_SNESVI *dmsnesvi1; 109 Vec finemarked,coarsemarked; 110 IS inactive; 111 VecScatter inject; 112 const PetscInt *index; 113 PetscInt n,k,cnt = 0,rstart,*coarseindex; 114 PetscScalar *marked; 115 116 PetscFunctionBegin; 117 ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject*)&isnes);CHKERRQ(ierr); 118 if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm1),PETSC_ERR_PLIB,"Composed VI data structure is missing"); 119 ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr); 120 121 /* get the original coarsen */ 122 ierr = (*dmsnesvi1->coarsen)(dm1,comm,dm2);CHKERRQ(ierr); 123 124 /* not sure why this extra reference is needed, but without the dm2 disappears too early */ 125 ierr = PetscObjectReference((PetscObject)*dm2);CHKERRQ(ierr); 126 127 /* need to set back global vectors in order to use the original injection */ 128 ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr); 129 130 dm1->ops->createglobalvector = dmsnesvi1->createglobalvector; 131 132 ierr = DMCreateGlobalVector(dm1,&finemarked);CHKERRQ(ierr); 133 ierr = DMCreateGlobalVector(*dm2,&coarsemarked);CHKERRQ(ierr); 134 135 /* 136 fill finemarked with locations of inactive points 137 */ 138 ierr = ISGetIndices(dmsnesvi1->inactive,&index);CHKERRQ(ierr); 139 ierr = ISGetLocalSize(dmsnesvi1->inactive,&n);CHKERRQ(ierr); 140 ierr = VecSet(finemarked,0.0);CHKERRQ(ierr); 141 for (k=0; k<n; k++) { 142 ierr = VecSetValue(finemarked,index[k],1.0,INSERT_VALUES);CHKERRQ(ierr); 143 } 144 ierr = VecAssemblyBegin(finemarked);CHKERRQ(ierr); 145 ierr = VecAssemblyEnd(finemarked);CHKERRQ(ierr); 146 147 ierr = DMCreateInjection(*dm2,dm1,&inject);CHKERRQ(ierr); 148 ierr = VecScatterBegin(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 149 ierr = VecScatterEnd(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 150 ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); 151 152 /* 153 create index set list of coarse inactive points from coarsemarked 154 */ 155 ierr = VecGetLocalSize(coarsemarked,&n);CHKERRQ(ierr); 156 ierr = VecGetOwnershipRange(coarsemarked,&rstart,NULL);CHKERRQ(ierr); 157 ierr = VecGetArray(coarsemarked,&marked);CHKERRQ(ierr); 158 for (k=0; k<n; k++) { 159 if (marked[k] != 0.0) cnt++; 160 } 161 ierr = PetscMalloc(cnt*sizeof(PetscInt),&coarseindex);CHKERRQ(ierr); 162 cnt = 0; 163 for (k=0; k<n; k++) { 164 if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart; 165 } 166 ierr = VecRestoreArray(coarsemarked,&marked);CHKERRQ(ierr); 167 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)coarsemarked),cnt,coarseindex,PETSC_OWN_POINTER,&inactive);CHKERRQ(ierr); 168 169 ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr); 170 171 dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI; 172 173 ierr = DMSetVI(*dm2,inactive);CHKERRQ(ierr); 174 175 ierr = VecDestroy(&finemarked);CHKERRQ(ierr); 176 ierr = VecDestroy(&coarsemarked);CHKERRQ(ierr); 177 ierr = ISDestroy(&inactive);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "DMDestroy_SNESVI" 183 PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi) 184 { 185 PetscErrorCode ierr; 186 187 PetscFunctionBegin; 188 /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */ 189 dmsnesvi->dm->ops->createinterpolation = dmsnesvi->createinterpolation; 190 dmsnesvi->dm->ops->coarsen = dmsnesvi->coarsen; 191 dmsnesvi->dm->ops->createglobalvector = dmsnesvi->createglobalvector; 192 /* need to clear out this vectors because some of them may not have a reference to the DM 193 but they are counted as having references to the DM in DMDestroy() */ 194 ierr = DMClearGlobalVectors(dmsnesvi->dm);CHKERRQ(ierr); 195 196 ierr = ISDestroy(&dmsnesvi->inactive);CHKERRQ(ierr); 197 ierr = PetscFree(dmsnesvi);CHKERRQ(ierr); 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "DMSetVI" 203 /* 204 DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to 205 be restricted to only those variables NOT associated with active constraints. 206 207 */ 208 PetscErrorCode DMSetVI(DM dm,IS inactive) 209 { 210 PetscErrorCode ierr; 211 PetscContainer isnes; 212 DM_SNESVI *dmsnesvi; 213 214 PetscFunctionBegin; 215 if (!dm) PetscFunctionReturn(0); 216 217 ierr = PetscObjectReference((PetscObject)inactive);CHKERRQ(ierr); 218 219 ierr = PetscObjectQuery((PetscObject)dm,"VI",(PetscObject*)&isnes);CHKERRQ(ierr); 220 if (!isnes) { 221 ierr = PetscContainerCreate(PetscObjectComm((PetscObject)dm),&isnes);CHKERRQ(ierr); 222 ierr = PetscContainerSetUserDestroy(isnes,(PetscErrorCode (*)(void*))DMDestroy_SNESVI);CHKERRQ(ierr); 223 ierr = PetscNew(DM_SNESVI,&dmsnesvi);CHKERRQ(ierr); 224 ierr = PetscContainerSetPointer(isnes,(void*)dmsnesvi);CHKERRQ(ierr); 225 ierr = PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)isnes);CHKERRQ(ierr); 226 ierr = PetscContainerDestroy(&isnes);CHKERRQ(ierr); 227 228 dmsnesvi->createinterpolation = dm->ops->createinterpolation; 229 dm->ops->createinterpolation = DMCreateInterpolation_SNESVI; 230 dmsnesvi->coarsen = dm->ops->coarsen; 231 dm->ops->coarsen = DMCoarsen_SNESVI; 232 dmsnesvi->createglobalvector = dm->ops->createglobalvector; 233 dm->ops->createglobalvector = DMCreateGlobalVector_SNESVI; 234 } else { 235 ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi);CHKERRQ(ierr); 236 ierr = ISDestroy(&dmsnesvi->inactive);CHKERRQ(ierr); 237 } 238 ierr = DMClearGlobalVectors(dm);CHKERRQ(ierr); 239 ierr = ISGetLocalSize(inactive,&dmsnesvi->n);CHKERRQ(ierr); 240 241 dmsnesvi->inactive = inactive; 242 dmsnesvi->dm = dm; 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "DMDestroyVI" 248 /* 249 DMDestroyVI - Frees the DM_SNESVI object contained in the DM 250 - also resets the function pointers in the DM for createinterpolation() etc to use the original DM 251 */ 252 PetscErrorCode DMDestroyVI(DM dm) 253 { 254 PetscErrorCode ierr; 255 256 PetscFunctionBegin; 257 if (!dm) PetscFunctionReturn(0); 258 ierr = PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)NULL);CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 /* --------------------------------------------------------------------------------------------------------*/ 263 264 265 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "SNESCreateIndexSets_VINEWTONRSLS" 269 PetscErrorCode SNESCreateIndexSets_VINEWTONRSLS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact) 270 { 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr); 275 ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 280 #undef __FUNCT__ 281 #define __FUNCT__ "SNESCreateSubVectors_VINEWTONRSLS" 282 PetscErrorCode SNESCreateSubVectors_VINEWTONRSLS(SNES snes,PetscInt n,Vec *newv) 283 { 284 PetscErrorCode ierr; 285 Vec v; 286 287 PetscFunctionBegin; 288 ierr = VecCreate(PetscObjectComm((PetscObject)snes),&v);CHKERRQ(ierr); 289 ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr); 290 ierr = VecSetFromOptions(v);CHKERRQ(ierr); 291 *newv = v; 292 PetscFunctionReturn(0); 293 } 294 295 /* Resets the snes PC and KSP when the active set sizes change */ 296 #undef __FUNCT__ 297 #define __FUNCT__ "SNESVIResetPCandKSP" 298 PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat) 299 { 300 PetscErrorCode ierr; 301 KSP snesksp; 302 303 PetscFunctionBegin; 304 ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr); 305 ierr = KSPReset(snesksp);CHKERRQ(ierr); 306 307 /* 308 KSP kspnew; 309 PC pcnew; 310 const MatSolverPackage stype; 311 312 313 ierr = KSPCreate(PetscObjectComm((PetscObject)snes),&kspnew);CHKERRQ(ierr); 314 kspnew->pc_side = snesksp->pc_side; 315 kspnew->rtol = snesksp->rtol; 316 kspnew->abstol = snesksp->abstol; 317 kspnew->max_it = snesksp->max_it; 318 ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr); 319 ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr); 320 ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr); 321 ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 322 ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr); 323 ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr); 324 ierr = KSPDestroy(&snesksp);CHKERRQ(ierr); 325 snes->ksp = kspnew; 326 ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr); 327 ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);*/ 328 PetscFunctionReturn(0); 329 } 330 331 /* Variational Inequality solver using reduce space method. No semismooth algorithm is 332 implemented in this algorithm. It basically identifies the active constraints and does 333 a linear solve on the other variables (those not associated with the active constraints). */ 334 #undef __FUNCT__ 335 #define __FUNCT__ "SNESSolve_VINEWTONRSLS" 336 PetscErrorCode SNESSolve_VINEWTONRSLS(SNES snes) 337 { 338 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*)snes->data; 339 PetscErrorCode ierr; 340 PetscInt maxits,i,lits; 341 PetscBool lssucceed; 342 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 343 PetscReal fnorm,gnorm,xnorm=0,ynorm; 344 Vec Y,X,F; 345 KSPConvergedReason kspreason; 346 347 PetscFunctionBegin; 348 snes->numFailures = 0; 349 snes->numLinearSolveFailures = 0; 350 snes->reason = SNES_CONVERGED_ITERATING; 351 352 maxits = snes->max_its; /* maximum number of iterations */ 353 X = snes->vec_sol; /* solution vector */ 354 F = snes->vec_func; /* residual vector */ 355 Y = snes->work[0]; /* work vectors */ 356 357 ierr = SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm);CHKERRQ(ierr); 358 ierr = SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 359 ierr = SNESLineSearchSetUp(snes->linesearch);CHKERRQ(ierr); 360 361 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 362 snes->iter = 0; 363 snes->norm = 0.0; 364 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 365 366 ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr); 367 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 368 if (snes->domainerror) { 369 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 370 PetscFunctionReturn(0); 371 } 372 ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr); 373 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 374 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 375 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 376 377 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 378 snes->norm = fnorm; 379 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 380 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 381 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 382 383 /* set parameter for default relative tolerance convergence test */ 384 snes->ttol = fnorm*snes->rtol; 385 /* test convergence */ 386 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 387 if (snes->reason) PetscFunctionReturn(0); 388 389 390 for (i=0; i<maxits; i++) { 391 392 IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 393 IS IS_redact; /* redundant active set */ 394 VecScatter scat_act,scat_inact; 395 PetscInt nis_act,nis_inact; 396 Vec Y_act,Y_inact,F_inact; 397 Mat jac_inact_inact,prejac_inact_inact; 398 PetscBool isequal; 399 400 /* Call general purpose update function */ 401 if (snes->ops->update) { 402 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 403 } 404 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 405 406 407 /* Create active and inactive index sets */ 408 409 /*original 410 ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr); 411 */ 412 ierr = SNESVIGetActiveSetIS(snes,X,F,&IS_act);CHKERRQ(ierr); 413 414 if (vi->checkredundancy) { 415 (*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);CHKERRQ(ierr); 416 if (IS_redact) { 417 ierr = ISSort(IS_redact);CHKERRQ(ierr); 418 ierr = ISComplement(IS_redact,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr); 419 ierr = ISDestroy(&IS_redact);CHKERRQ(ierr); 420 } else { 421 ierr = ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr); 422 } 423 } else { 424 ierr = ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr); 425 } 426 427 428 /* Create inactive set submatrix */ 429 ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 430 431 if (0) { /* Dead code (temporary developer hack) */ 432 IS keptrows; 433 ierr = MatFindNonzeroRows(jac_inact_inact,&keptrows);CHKERRQ(ierr); 434 if (keptrows) { 435 PetscInt cnt,*nrows,k; 436 const PetscInt *krows,*inact; 437 PetscInt rstart=jac_inact_inact->rmap->rstart; 438 439 ierr = MatDestroy(&jac_inact_inact);CHKERRQ(ierr); 440 ierr = ISDestroy(&IS_act);CHKERRQ(ierr); 441 442 ierr = ISGetLocalSize(keptrows,&cnt);CHKERRQ(ierr); 443 ierr = ISGetIndices(keptrows,&krows);CHKERRQ(ierr); 444 ierr = ISGetIndices(IS_inact,&inact);CHKERRQ(ierr); 445 ierr = PetscMalloc(cnt*sizeof(PetscInt),&nrows);CHKERRQ(ierr); 446 for (k=0; k<cnt; k++) nrows[k] = inact[krows[k]-rstart]; 447 ierr = ISRestoreIndices(keptrows,&krows);CHKERRQ(ierr); 448 ierr = ISRestoreIndices(IS_inact,&inact);CHKERRQ(ierr); 449 ierr = ISDestroy(&keptrows);CHKERRQ(ierr); 450 ierr = ISDestroy(&IS_inact);CHKERRQ(ierr); 451 452 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),cnt,nrows,PETSC_OWN_POINTER,&IS_inact);CHKERRQ(ierr); 453 ierr = ISComplement(IS_inact,F->map->rstart,F->map->rend,&IS_act);CHKERRQ(ierr); 454 ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 455 } 456 } 457 ierr = DMSetVI(snes->dm,IS_inact);CHKERRQ(ierr); 458 /* remove later */ 459 460 /* 461 ierr = VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm));CHKERRQ(ierr); 462 ierr = VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm));CHKERRQ(ierr); 463 ierr = VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X)));CHKERRQ(ierr); 464 ierr = VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F)));CHKERRQ(ierr); 465 ierr = ISView(IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)IS_inact)));CHKERRQ(ierr); 466 */ 467 468 /* Get sizes of active and inactive sets */ 469 ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 470 ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 471 472 /* Create active and inactive set vectors */ 473 ierr = SNESCreateSubVectors_VINEWTONRSLS(snes,nis_inact,&F_inact);CHKERRQ(ierr); 474 ierr = SNESCreateSubVectors_VINEWTONRSLS(snes,nis_act,&Y_act);CHKERRQ(ierr); 475 ierr = SNESCreateSubVectors_VINEWTONRSLS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 476 477 /* Create scatter contexts */ 478 ierr = VecScatterCreate(Y,IS_act,Y_act,NULL,&scat_act);CHKERRQ(ierr); 479 ierr = VecScatterCreate(Y,IS_inact,Y_inact,NULL,&scat_inact);CHKERRQ(ierr); 480 481 /* Do a vec scatter to active and inactive set vectors */ 482 ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 483 ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 484 485 ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 486 ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 487 488 ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 489 ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 490 491 /* Active set direction = 0 */ 492 ierr = VecSet(Y_act,0);CHKERRQ(ierr); 493 if (snes->jacobian != snes->jacobian_pre) { 494 ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 495 } else prejac_inact_inact = jac_inact_inact; 496 497 ierr = ISEqual(vi->IS_inact_prev,IS_inact,&isequal);CHKERRQ(ierr); 498 if (!isequal) { 499 ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr); 500 flg = DIFFERENT_NONZERO_PATTERN; 501 } 502 503 /* ierr = ISView(IS_inact,0);CHKERRQ(ierr); */ 504 /* ierr = ISView(IS_act,0);CHKERRQ(ierr);*/ 505 /* ierr = MatView(snes->jacobian_pre,0); */ 506 507 508 509 ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 510 ierr = KSPSetUp(snes->ksp);CHKERRQ(ierr); 511 { 512 PC pc; 513 PetscBool flg; 514 ierr = KSPGetPC(snes->ksp,&pc);CHKERRQ(ierr); 515 ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);CHKERRQ(ierr); 516 if (flg) { 517 KSP *subksps; 518 ierr = PCFieldSplitGetSubKSP(pc,NULL,&subksps);CHKERRQ(ierr); 519 ierr = KSPGetPC(subksps[0],&pc);CHKERRQ(ierr); 520 ierr = PetscFree(subksps);CHKERRQ(ierr); 521 ierr = PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&flg);CHKERRQ(ierr); 522 if (flg) { 523 PetscInt n,N = 101*101,j,cnts[3] = {0,0,0}; 524 const PetscInt *ii; 525 526 ierr = ISGetSize(IS_inact,&n);CHKERRQ(ierr); 527 ierr = ISGetIndices(IS_inact,&ii);CHKERRQ(ierr); 528 for (j=0; j<n; j++) { 529 if (ii[j] < N) cnts[0]++; 530 else if (ii[j] < 2*N) cnts[1]++; 531 else if (ii[j] < 3*N) cnts[2]++; 532 } 533 ierr = ISRestoreIndices(IS_inact,&ii);CHKERRQ(ierr); 534 535 ierr = PCBJacobiSetTotalBlocks(pc,3,cnts);CHKERRQ(ierr); 536 } 537 } 538 } 539 540 ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr); 541 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 542 if (kspreason < 0) { 543 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 544 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 545 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 546 break; 547 } 548 } 549 550 ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 551 ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 552 ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 553 ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 554 555 ierr = VecDestroy(&F_inact);CHKERRQ(ierr); 556 ierr = VecDestroy(&Y_act);CHKERRQ(ierr); 557 ierr = VecDestroy(&Y_inact);CHKERRQ(ierr); 558 ierr = VecScatterDestroy(&scat_act);CHKERRQ(ierr); 559 ierr = VecScatterDestroy(&scat_inact);CHKERRQ(ierr); 560 ierr = ISDestroy(&IS_act);CHKERRQ(ierr); 561 if (!isequal) { 562 ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr); 563 ierr = ISDuplicate(IS_inact,&vi->IS_inact_prev);CHKERRQ(ierr); 564 } 565 ierr = ISDestroy(&IS_inact);CHKERRQ(ierr); 566 ierr = MatDestroy(&jac_inact_inact);CHKERRQ(ierr); 567 if (snes->jacobian != snes->jacobian_pre) { 568 ierr = MatDestroy(&prejac_inact_inact);CHKERRQ(ierr); 569 } 570 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 571 snes->linear_its += lits; 572 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 573 /* 574 if (snes->ops->precheck) { 575 PetscBool changed_y = PETSC_FALSE; 576 ierr = (*snes->ops->precheck)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr); 577 } 578 579 if (PetscLogPrintInfo) { 580 ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 581 } 582 */ 583 /* Compute a (scaled) negative update in the line search routine: 584 Y <- X - lambda*Y 585 and evaluate G = function(Y) (depends on the line search). 586 */ 587 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 588 ynorm = 1; gnorm = fnorm; 589 ierr = SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y);CHKERRQ(ierr); 590 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 591 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 592 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 593 if (snes->domainerror) { 594 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 595 ierr = DMDestroyVI(snes->dm);CHKERRQ(ierr); 596 PetscFunctionReturn(0); 597 } 598 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 599 600 if (!lssucceed) { 601 if (++snes->numFailures >= snes->maxFailures) { 602 PetscBool ismin; 603 snes->reason = SNES_DIVERGED_LINE_SEARCH; 604 ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,F,X,gnorm,&ismin);CHKERRQ(ierr); 605 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 606 break; 607 } 608 } 609 /* Update function and solution vectors */ 610 fnorm = gnorm; 611 /* Monitor convergence */ 612 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 613 snes->iter = i+1; 614 snes->norm = fnorm; 615 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 616 ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 617 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 618 /* Test for convergence, xnorm = || X || */ 619 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 620 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 621 if (snes->reason) break; 622 } 623 ierr = DMDestroyVI(snes->dm);CHKERRQ(ierr); 624 if (i == maxits) { 625 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 626 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 627 } 628 PetscFunctionReturn(0); 629 } 630 631 #undef __FUNCT__ 632 #define __FUNCT__ "SNESVISetRedundancyCheck" 633 PetscErrorCode SNESVISetRedundancyCheck(SNES snes,PetscErrorCode (*func)(SNES,IS,IS*,void*),void *ctx) 634 { 635 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*)snes->data; 636 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 639 vi->checkredundancy = func; 640 vi->ctxP = ctx; 641 PetscFunctionReturn(0); 642 } 643 644 #if defined(PETSC_HAVE_MATLAB_ENGINE) 645 #include <engine.h> 646 #include <mex.h> 647 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext; 648 649 #undef __FUNCT__ 650 #define __FUNCT__ "SNESVIRedundancyCheck_Matlab" 651 PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes,IS is_act,IS *is_redact,void *ctx) 652 { 653 PetscErrorCode ierr; 654 SNESMatlabContext *sctx = (SNESMatlabContext*)ctx; 655 int nlhs = 1, nrhs = 5; 656 mxArray *plhs[1], *prhs[5]; 657 long long int l1 = 0, l2 = 0, ls = 0; 658 PetscInt *indices=NULL; 659 660 PetscFunctionBegin; 661 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 662 PetscValidHeaderSpecific(is_act,IS_CLASSID,2); 663 PetscValidPointer(is_redact,3); 664 PetscCheckSameComm(snes,1,is_act,2); 665 666 /* Create IS for reduced active set of size 0, its size and indices will 667 bet set by the Matlab function */ 668 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),0,indices,PETSC_OWN_POINTER,is_redact);CHKERRQ(ierr); 669 /* call Matlab function in ctx */ 670 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 671 ierr = PetscMemcpy(&l1,&is_act,sizeof(is_act));CHKERRQ(ierr); 672 ierr = PetscMemcpy(&l2,is_redact,sizeof(is_act));CHKERRQ(ierr); 673 prhs[0] = mxCreateDoubleScalar((double)ls); 674 prhs[1] = mxCreateDoubleScalar((double)l1); 675 prhs[2] = mxCreateDoubleScalar((double)l2); 676 prhs[3] = mxCreateString(sctx->funcname); 677 prhs[4] = sctx->ctx; 678 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESVIRedundancyCheckInternal");CHKERRQ(ierr); 679 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 680 mxDestroyArray(prhs[0]); 681 mxDestroyArray(prhs[1]); 682 mxDestroyArray(prhs[2]); 683 mxDestroyArray(prhs[3]); 684 mxDestroyArray(plhs[0]); 685 PetscFunctionReturn(0); 686 } 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "SNESVISetRedundancyCheckMatlab" 690 PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes,const char *func,mxArray *ctx) 691 { 692 PetscErrorCode ierr; 693 SNESMatlabContext *sctx; 694 695 PetscFunctionBegin; 696 /* currently sctx is memory bleed */ 697 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 698 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 699 sctx->ctx = mxDuplicateArray(ctx); 700 ierr = SNESVISetRedundancyCheck(snes,SNESVIRedundancyCheck_Matlab,sctx);CHKERRQ(ierr); 701 PetscFunctionReturn(0); 702 } 703 704 #endif 705 706 /* -------------------------------------------------------------------------- */ 707 /* 708 SNESSetUp_VINEWTONRSLS - Sets up the internal data structures for the later use 709 of the SNESVI nonlinear solver. 710 711 Input Parameter: 712 . snes - the SNES context 713 . x - the solution vector 714 715 Application Interface Routine: SNESSetUp() 716 717 Notes: 718 For basic use of the SNES solvers, the user need not explicitly call 719 SNESSetUp(), since these actions will automatically occur during 720 the call to SNESSolve(). 721 */ 722 #undef __FUNCT__ 723 #define __FUNCT__ "SNESSetUp_VINEWTONRSLS" 724 PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes) 725 { 726 PetscErrorCode ierr; 727 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*) snes->data; 728 PetscInt *indices; 729 PetscInt i,n,rstart,rend; 730 SNESLineSearch linesearch; 731 732 PetscFunctionBegin; 733 ierr = SNESSetUp_VI(snes);CHKERRQ(ierr); 734 735 /* Set up previous active index set for the first snes solve 736 vi->IS_inact_prev = 0,1,2,....N */ 737 738 ierr = VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);CHKERRQ(ierr); 739 ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 740 ierr = PetscMalloc(n*sizeof(PetscInt),&indices);CHKERRQ(ierr); 741 for (i=0; i < n; i++) indices[i] = rstart + i; 742 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);CHKERRQ(ierr); 743 744 /* set the line search functions */ 745 if (!snes->linesearch) { 746 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 747 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 748 } 749 PetscFunctionReturn(0); 750 } 751 /* -------------------------------------------------------------------------- */ 752 #undef __FUNCT__ 753 #define __FUNCT__ "SNESReset_VINEWTONRSLS" 754 PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes) 755 { 756 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*) snes->data; 757 PetscErrorCode ierr; 758 759 PetscFunctionBegin; 760 ierr = SNESReset_VI(snes);CHKERRQ(ierr); 761 ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr); 762 PetscFunctionReturn(0); 763 } 764 765 /* -------------------------------------------------------------------------- */ 766 /*MC 767 SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method 768 769 Options Database: 770 + -snes_vi_type <ss,rs,rsaug> a semi-smooth solver, a reduced space active set method, and a reduced space active set method that does not eliminate the active constraints from the Jacobian instead augments the Jacobian with additional variables that enforce the constraints 771 - -snes_vi_monitor - prints the number of active constraints at each iteration. 772 773 Level: beginner 774 775 References: 776 - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large-Scale 777 Applications, Optimization Methods and Software, 21 (2006). 778 779 .seealso: SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVINEWTONRSLS, SNESVINEWTONSSLS, SNESNEWTONTR, SNESLineSearchSet(), 780 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 781 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 782 783 M*/ 784 EXTERN_C_BEGIN 785 #undef __FUNCT__ 786 #define __FUNCT__ "SNESCreate_VINEWTONRSLS" 787 PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes) 788 { 789 PetscErrorCode ierr; 790 SNES_VINEWTONRSLS *vi; 791 792 PetscFunctionBegin; 793 snes->ops->reset = SNESReset_VINEWTONRSLS; 794 snes->ops->setup = SNESSetUp_VINEWTONRSLS; 795 snes->ops->solve = SNESSolve_VINEWTONRSLS; 796 snes->ops->destroy = SNESDestroy_VI; 797 snes->ops->setfromoptions = SNESSetFromOptions_VI; 798 snes->ops->view = NULL; 799 snes->ops->converged = SNESDefaultConverged_VI; 800 801 snes->usesksp = PETSC_TRUE; 802 snes->usespc = PETSC_FALSE; 803 804 ierr = PetscNewLog(snes,SNES_VINEWTONRSLS,&vi);CHKERRQ(ierr); 805 snes->data = (void*)vi; 806 vi->checkredundancy = NULL; 807 808 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetVariableBounds_C","SNESVISetVariableBounds_VI",SNESVISetVariableBounds_VI);CHKERRQ(ierr); 809 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetComputeVariableBounds_C","SNESVISetComputeVariableBounds_VI",SNESVISetComputeVariableBounds_VI);CHKERRQ(ierr); 810 PetscFunctionReturn(0); 811 } 812 EXTERN_C_END 813 814