1 #include <petscdmda.h> /*I "petscdmda.h" I*/ 2 #include <petsc/private/dmimpl.h> 3 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 4 5 /* This structure holds the user-provided DMDA callbacks */ 6 typedef struct { 7 PetscErrorCode (*residuallocal)(DMDALocalInfo*,void*,void*,void*); 8 PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,void*); 9 PetscErrorCode (*objectivelocal)(DMDALocalInfo*,void*,PetscReal*,void*); 10 11 PetscErrorCode (*residuallocalvec)(DMDALocalInfo*,Vec,Vec,void*); 12 PetscErrorCode (*jacobianlocalvec)(DMDALocalInfo*,Vec,Mat,Mat,void*); 13 PetscErrorCode (*objectivelocalvec)(DMDALocalInfo*,Vec,PetscReal*,void*); 14 void *residuallocalctx; 15 void *jacobianlocalctx; 16 void *objectivelocalctx; 17 InsertMode residuallocalimode; 18 19 /* For Picard iteration defined locally */ 20 PetscErrorCode (*rhsplocal)(DMDALocalInfo*,void*,void*,void*); 21 PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,void*); 22 void *picardlocalctx; 23 } DMSNES_DA; 24 25 static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) 26 { 27 PetscFunctionBegin; 28 PetscCall(PetscFree(sdm->data)); 29 PetscFunctionReturn(0); 30 } 31 32 static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm) 33 { 34 PetscFunctionBegin; 35 PetscCall(PetscNewLog(sdm,(DMSNES_DA**)&sdm->data)); 36 if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA))); 37 PetscFunctionReturn(0); 38 } 39 40 static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA **dmdasnes) 41 { 42 PetscFunctionBegin; 43 *dmdasnes = NULL; 44 if (!sdm->data) { 45 PetscCall(PetscNewLog(dm,(DMSNES_DA**)&sdm->data)); 46 sdm->ops->destroy = DMSNESDestroy_DMDA; 47 sdm->ops->duplicate = DMSNESDuplicate_DMDA; 48 } 49 *dmdasnes = (DMSNES_DA*)sdm->data; 50 PetscFunctionReturn(0); 51 } 52 53 static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 54 { 55 DM dm; 56 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 57 DMDALocalInfo info; 58 Vec Xloc; 59 void *x,*f; 60 61 PetscFunctionBegin; 62 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 63 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 64 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 65 PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 66 PetscCall(SNESGetDM(snes,&dm)); 67 PetscCall(DMGetLocalVector(dm,&Xloc)); 68 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 69 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 70 PetscCall(DMDAGetLocalInfo(dm,&info)); 71 switch (dmdasnes->residuallocalimode) { 72 case INSERT_VALUES: { 73 PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0)); 74 if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocalvec)(&info,Xloc,F,dmdasnes->residuallocalctx)); 75 else { 76 PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 77 PetscCall(DMDAVecGetArray(dm,F,&f)); 78 PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx)); 79 PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 80 PetscCall(DMDAVecRestoreArray(dm,F,&f)); 81 } 82 PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0)); 83 } break; 84 case ADD_VALUES: { 85 Vec Floc; 86 PetscCall(DMGetLocalVector(dm,&Floc)); 87 PetscCall(VecZeroEntries(Floc)); 88 PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0)); 89 if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocalvec)(&info,Xloc,Floc,dmdasnes->residuallocalctx)); 90 else { 91 PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 92 PetscCall(DMDAVecGetArray(dm,Floc,&f)); 93 PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx)); 94 PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 95 PetscCall(DMDAVecRestoreArray(dm,Floc,&f)); 96 } 97 PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0)); 98 PetscCall(VecZeroEntries(F)); 99 PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 100 PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 101 PetscCall(DMRestoreLocalVector(dm,&Floc)); 102 } break; 103 default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 104 } 105 PetscCall(DMRestoreLocalVector(dm,&Xloc)); 106 if (snes->domainerror) PetscCall(VecSetInf(F)); 107 PetscFunctionReturn(0); 108 } 109 110 static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx) 111 { 112 DM dm; 113 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 114 DMDALocalInfo info; 115 Vec Xloc; 116 void *x; 117 118 PetscFunctionBegin; 119 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 120 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 121 PetscValidRealPointer(ob,3); 122 PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 123 PetscCall(SNESGetDM(snes,&dm)); 124 PetscCall(DMGetLocalVector(dm,&Xloc)); 125 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 126 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 127 PetscCall(DMDAGetLocalInfo(dm,&info)); 128 if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective",(*dmdasnes->objectivelocalvec)(&info,Xloc,ob,dmdasnes->objectivelocalctx)); 129 else { 130 PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 131 PetscCallBack("SNES DMDA local callback objective",(*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx)); 132 PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 133 } 134 PetscCall(DMRestoreLocalVector(dm,&Xloc)); 135 PetscFunctionReturn(0); 136 } 137 138 /* Routine is called by example, hence must be labeled PETSC_EXTERN */ 139 PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 140 { 141 DM dm; 142 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 143 DMDALocalInfo info; 144 Vec Xloc; 145 void *x; 146 147 PetscFunctionBegin; 148 PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 149 PetscCall(SNESGetDM(snes,&dm)); 150 151 if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) { 152 PetscCall(DMGetLocalVector(dm,&Xloc)); 153 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 154 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 155 PetscCall(DMDAGetLocalInfo(dm,&info)); 156 if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian",(*dmdasnes->jacobianlocalvec)(&info,Xloc,A,B,dmdasnes->jacobianlocalctx)); 157 else { 158 PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 159 PetscCallBack("SNES DMDA local callback Jacobian",(*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx)); 160 PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 161 } 162 PetscCall(DMRestoreLocalVector(dm,&Xloc)); 163 } else { 164 MatFDColoring fdcoloring; 165 PetscCall(PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring)); 166 if (!fdcoloring) { 167 ISColoring coloring; 168 169 PetscCall(DMCreateColoring(dm,dm->coloringtype,&coloring)); 170 PetscCall(MatFDColoringCreate(B,coloring,&fdcoloring)); 171 switch (dm->coloringtype) { 172 case IS_COLORING_GLOBAL: 173 PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes)); 174 break; 175 default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 176 } 177 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix)); 178 PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 179 PetscCall(MatFDColoringSetUp(B,coloring,fdcoloring)); 180 PetscCall(ISColoringDestroy(&coloring)); 181 PetscCall(PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring)); 182 PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 183 184 /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 185 * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 186 * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 187 * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 188 * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 189 */ 190 PetscCall(PetscObjectDereference((PetscObject)dm)); 191 } 192 PetscCall(MatFDColoringApply(B,fdcoloring,X,snes)); 193 } 194 /* This will be redundant if the user called both, but it's too common to forget. */ 195 if (A != B) { 196 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 197 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 198 } 199 PetscFunctionReturn(0); 200 } 201 202 /*@C 203 DMDASNESSetFunctionLocal - set a local residual evaluation function 204 205 Logically Collective 206 207 Input Parameters: 208 + dm - DM to associate callback with 209 . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 210 . func - local residual evaluation 211 - ctx - optional context for local residual evaluation 212 213 Calling sequence: 214 For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx), 215 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 216 . x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x) 217 . f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f) 218 - ctx - optional context passed above 219 220 Level: beginner 221 222 .seealso: `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 223 @*/ 224 PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 225 { 226 DMSNES sdm; 227 DMSNES_DA *dmdasnes; 228 229 PetscFunctionBegin; 230 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 231 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 232 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 233 234 dmdasnes->residuallocalimode = imode; 235 dmdasnes->residuallocal = func; 236 dmdasnes->residuallocalctx = ctx; 237 238 PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes)); 239 if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 240 PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 241 } 242 PetscFunctionReturn(0); 243 } 244 245 /*@C 246 DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector 247 248 Logically Collective 249 250 Input Parameters: 251 + dm - DM to associate callback with 252 . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 253 . func - local residual evaluation 254 - ctx - optional context for local residual evaluation 255 256 Calling sequence: 257 For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x, Vec f, void *ctx), 258 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 259 . x - state vector at which to evaluate residual 260 . f - residual vector 261 - ctx - optional context passed above 262 263 Level: beginner 264 265 .seealso: `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 266 @*/ 267 PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Vec,void*),void *ctx) 268 { 269 DMSNES sdm; 270 DMSNES_DA *dmdasnes; 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 274 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 275 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 276 277 dmdasnes->residuallocalimode = imode; 278 dmdasnes->residuallocalvec = func; 279 dmdasnes->residuallocalctx = ctx; 280 281 PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes)); 282 if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 283 PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 284 } 285 PetscFunctionReturn(0); 286 } 287 288 /*@C 289 DMDASNESSetJacobianLocal - set a local Jacobian evaluation function 290 291 Logically Collective 292 293 Input Parameters: 294 + dm - DM to associate callback with 295 . func - local Jacobian evaluation 296 - ctx - optional context for local Jacobian evaluation 297 298 Calling sequence: 299 For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx), 300 + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 301 . x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x) 302 . J - Mat object for the Jacobian 303 . M - Mat object for the Jacobian preconditioner matrix 304 - ctx - optional context passed above 305 306 Level: beginner 307 308 .seealso: `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 309 @*/ 310 PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 311 { 312 DMSNES sdm; 313 DMSNES_DA *dmdasnes; 314 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 317 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 318 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 319 320 dmdasnes->jacobianlocal = func; 321 dmdasnes->jacobianlocalctx = ctx; 322 323 PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 324 PetscFunctionReturn(0); 325 } 326 327 /*@C 328 DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector 329 330 Logically Collective 331 332 Input Parameters: 333 + dm - DM to associate callback with 334 . func - local Jacobian evaluation 335 - ctx - optional context for local Jacobian evaluation 336 337 Calling sequence: 338 For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,Mat J,Mat M,void *ctx), 339 + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 340 . x - state vector at which to evaluate Jacobian 341 . J - Mat object for the Jacobian 342 . M - Mat object for the Jacobian preconditioner matrix 343 - ctx - optional context passed above 344 345 Level: beginner 346 347 .seealso: `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 348 @*/ 349 PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Mat,Mat,void*),void *ctx) 350 { 351 DMSNES sdm; 352 DMSNES_DA *dmdasnes; 353 354 PetscFunctionBegin; 355 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 356 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 357 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 358 359 dmdasnes->jacobianlocalvec = func; 360 dmdasnes->jacobianlocalctx = ctx; 361 362 PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 363 PetscFunctionReturn(0); 364 } 365 366 /*@C 367 DMDASNESSetObjectiveLocal - set a local residual evaluation function 368 369 Logically Collective 370 371 Input Parameters: 372 + dm - DM to associate callback with 373 . func - local objective evaluation 374 - ctx - optional context for local residual evaluation 375 376 Calling sequence for func: 377 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 378 . x - dimensional pointer to state at which to evaluate residual 379 . ob - eventual objective value 380 - ctx - optional context passed above 381 382 Level: beginner 383 384 .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 385 @*/ 386 PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx) 387 { 388 DMSNES sdm; 389 DMSNES_DA *dmdasnes; 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 393 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 394 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 395 396 dmdasnes->objectivelocal = func; 397 dmdasnes->objectivelocalctx = ctx; 398 399 PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes)); 400 PetscFunctionReturn(0); 401 } 402 403 /*@C 404 DMDASNESSetObjectiveLocal - set a local residual evaluation function that operates on a local vector 405 406 Logically Collective 407 408 Input Parameters: 409 + dm - DM to associate callback with 410 . func - local objective evaluation 411 - ctx - optional context for local residual evaluation 412 413 Calling sequence 414 For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,PetscReal *ob,void *ctx); 415 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 416 . x - state vector at which to evaluate residual 417 . ob - eventual objective value 418 - ctx - optional context passed above 419 420 Level: beginner 421 422 .seealso: `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 423 @*/ 424 PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm,DMDASNESObjectiveVec func,void *ctx) 425 { 426 DMSNES sdm; 427 DMSNES_DA *dmdasnes; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 431 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 432 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 433 434 dmdasnes->objectivelocalvec = func; 435 dmdasnes->objectivelocalctx = ctx; 436 437 PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes)); 438 PetscFunctionReturn(0); 439 } 440 441 static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx) 442 { 443 DM dm; 444 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 445 DMDALocalInfo info; 446 Vec Xloc; 447 void *x,*f; 448 449 PetscFunctionBegin; 450 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 451 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 452 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 453 PetscCheck(dmdasnes->rhsplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 454 PetscCall(SNESGetDM(snes,&dm)); 455 PetscCall(DMGetLocalVector(dm,&Xloc)); 456 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 457 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 458 PetscCall(DMDAGetLocalInfo(dm,&info)); 459 PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 460 switch (dmdasnes->residuallocalimode) { 461 case INSERT_VALUES: { 462 PetscCall(DMDAVecGetArray(dm,F,&f)); 463 PetscCallBack("SNES Picard DMDA local callback function",(*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx)); 464 PetscCall(DMDAVecRestoreArray(dm,F,&f)); 465 } break; 466 case ADD_VALUES: { 467 Vec Floc; 468 PetscCall(DMGetLocalVector(dm,&Floc)); 469 PetscCall(VecZeroEntries(Floc)); 470 PetscCall(DMDAVecGetArray(dm,Floc,&f)); 471 PetscCallBack("SNES Picard DMDA local callback function",(*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx)); 472 PetscCall(DMDAVecRestoreArray(dm,Floc,&f)); 473 PetscCall(VecZeroEntries(F)); 474 PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 475 PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 476 PetscCall(DMRestoreLocalVector(dm,&Floc)); 477 } break; 478 default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 479 } 480 PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 481 PetscCall(DMRestoreLocalVector(dm,&Xloc)); 482 PetscFunctionReturn(0); 483 } 484 485 static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 486 { 487 DM dm; 488 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 489 DMDALocalInfo info; 490 Vec Xloc; 491 void *x; 492 493 PetscFunctionBegin; 494 PetscCheck(dmdasnes->jacobianplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 495 PetscCall(SNESGetDM(snes,&dm)); 496 497 PetscCall(DMGetLocalVector(dm,&Xloc)); 498 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 499 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 500 PetscCall(DMDAGetLocalInfo(dm,&info)); 501 PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 502 PetscCallBack("SNES Picard DMDA local callback Jacobian",(*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx)); 503 PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 504 PetscCall(DMRestoreLocalVector(dm,&Xloc)); 505 if (A != B) { 506 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 507 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 508 } 509 PetscFunctionReturn(0); 510 } 511 512 /*@C 513 DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration 514 515 Logically Collective 516 517 Input Parameters: 518 + dm - DM to associate callback with 519 . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 520 . func - local residual evaluation 521 - ctx - optional context for local residual evaluation 522 523 Calling sequence for func: 524 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 525 . x - dimensional pointer to state at which to evaluate residual 526 . f - dimensional pointer to residual, write the residual here 527 - ctx - optional context passed above 528 529 Notes: 530 The user must use 531 PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user)); 532 in their code before calling this routine. 533 534 Level: beginner 535 536 .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 537 @*/ 538 PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*), 539 PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 540 { 541 DMSNES sdm; 542 DMSNES_DA *dmdasnes; 543 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 546 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 547 PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 548 549 dmdasnes->residuallocalimode = imode; 550 dmdasnes->rhsplocal = func; 551 dmdasnes->jacobianplocal = jac; 552 dmdasnes->picardlocalctx = ctx; 553 554 PetscCall(DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes)); 555 PetscCall(DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes)); 556 PetscFunctionReturn(0); 557 } 558