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