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