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