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 if (snes->domainerror) { 113 ierr = VecSetInf(F);CHKERRQ(ierr); 114 } 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "SNESComputeObjective_DMDA" 120 static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx) 121 { 122 PetscErrorCode ierr; 123 DM dm; 124 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 125 DMDALocalInfo info; 126 Vec Xloc; 127 void *x; 128 129 PetscFunctionBegin; 130 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 131 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 132 PetscValidPointer(ob,3); 133 if (!dmdasnes->objectivelocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 134 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 135 ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 136 ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 137 ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 138 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 139 ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 140 CHKMEMQ; 141 ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr); 142 CHKMEMQ; 143 ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 144 ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 148 149 #undef __FUNCT__ 150 #define __FUNCT__ "SNESComputeJacobian_DMDA" 151 PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 152 { 153 PetscErrorCode ierr; 154 DM dm; 155 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 156 DMDALocalInfo info; 157 Vec Xloc; 158 void *x; 159 160 PetscFunctionBegin; 161 if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 162 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 163 164 if (dmdasnes->jacobianlocal) { 165 ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 166 ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 167 ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 168 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 169 ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 170 CHKMEMQ; 171 ierr = (*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 172 CHKMEMQ; 173 ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 174 ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 175 } else { 176 MatFDColoring fdcoloring; 177 ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 178 if (!fdcoloring) { 179 ISColoring coloring; 180 181 ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 182 ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 183 switch (dm->coloringtype) { 184 case IS_COLORING_GLOBAL: 185 ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 186 break; 187 default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 188 } 189 ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 190 ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 191 ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 192 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 193 ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 194 ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 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 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 203 } 204 ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 205 } 206 /* This will be redundant if the user called both, but it's too common to forget. */ 207 if (A != B) { 208 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 209 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 210 } 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "DMDASNESSetFunctionLocal" 216 /*@C 217 DMDASNESSetFunctionLocal - set a local residual evaluation function 218 219 Logically Collective 220 221 Input Arguments: 222 + dm - DM to associate callback with 223 . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 224 . func - local residual evaluation 225 - ctx - optional context for local residual evaluation 226 227 Calling sequence: 228 For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx), 229 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 230 . x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x) 231 . f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f) 232 - ctx - optional context passed above 233 234 Level: beginner 235 236 .seealso: DMDASNESSetJacobianLocal(), DMSNESSetFunction(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 237 @*/ 238 PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 239 { 240 PetscErrorCode ierr; 241 DMSNES sdm; 242 DMSNES_DA *dmdasnes; 243 244 PetscFunctionBegin; 245 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 246 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 247 ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 248 249 dmdasnes->residuallocalimode = imode; 250 dmdasnes->residuallocal = func; 251 dmdasnes->residuallocalctx = ctx; 252 253 ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 254 if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 255 ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 256 } 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "DMDASNESSetJacobianLocal" 262 /*@C 263 DMDASNESSetJacobianLocal - set a local Jacobian evaluation function 264 265 Logically Collective 266 267 Input Arguments: 268 + dm - DM to associate callback with 269 . func - local Jacobian evaluation 270 - ctx - optional context for local Jacobian evaluation 271 272 Calling sequence: 273 For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx), 274 + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 275 . x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x) 276 . J - Mat object for the Jacobian 277 . M - Mat object for the Jacobian preconditioner matrix 278 - ctx - optional context passed above 279 280 Level: beginner 281 282 .seealso: DMDASNESSetFunctionLocal(), DMSNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 283 @*/ 284 PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 285 { 286 PetscErrorCode ierr; 287 DMSNES sdm; 288 DMSNES_DA *dmdasnes; 289 290 PetscFunctionBegin; 291 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 292 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 293 ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 294 295 dmdasnes->jacobianlocal = func; 296 dmdasnes->jacobianlocalctx = ctx; 297 298 ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 299 PetscFunctionReturn(0); 300 } 301 302 303 #undef __FUNCT__ 304 #define __FUNCT__ "DMDASNESSetObjectiveLocal" 305 /*@C 306 DMDASNESSetObjectiveLocal - set a local residual evaluation function 307 308 Logically Collective 309 310 Input Arguments: 311 + dm - DM to associate callback with 312 . func - local objective evaluation 313 - ctx - optional context for local residual evaluation 314 315 Calling sequence for func: 316 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 317 . x - dimensional pointer to state at which to evaluate residual 318 . ob - eventual objective value 319 - ctx - optional context passed above 320 321 Level: beginner 322 323 .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 324 @*/ 325 PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx) 326 { 327 PetscErrorCode ierr; 328 DMSNES sdm; 329 DMSNES_DA *dmdasnes; 330 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 333 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 334 ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 335 336 dmdasnes->objectivelocal = func; 337 dmdasnes->objectivelocalctx = ctx; 338 339 ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr); 340 PetscFunctionReturn(0); 341 } 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "SNESComputePicard_DMDA" 345 static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx) 346 { 347 PetscErrorCode ierr; 348 DM dm; 349 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 350 DMDALocalInfo info; 351 Vec Xloc; 352 void *x,*f; 353 354 PetscFunctionBegin; 355 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 356 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 357 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 358 if (!dmdasnes->rhsplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 359 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 360 ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 361 ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 362 ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 363 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 364 ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 365 switch (dmdasnes->residuallocalimode) { 366 case INSERT_VALUES: { 367 ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 368 CHKMEMQ; 369 ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 370 CHKMEMQ; 371 ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 372 } break; 373 case ADD_VALUES: { 374 Vec Floc; 375 ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 376 ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 377 ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 378 CHKMEMQ; 379 ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr); 380 CHKMEMQ; 381 ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 382 ierr = VecZeroEntries(F);CHKERRQ(ierr); 383 ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 384 ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 385 ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 386 } break; 387 default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 388 } 389 ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 390 ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 391 PetscFunctionReturn(0); 392 } 393 394 #undef __FUNCT__ 395 #define __FUNCT__ "SNESComputePicardJacobian_DMDA" 396 static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 397 { 398 PetscErrorCode ierr; 399 DM dm; 400 DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 401 DMDALocalInfo info; 402 Vec Xloc; 403 void *x; 404 405 PetscFunctionBegin; 406 if (!dmdasnes->jacobianplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 407 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 408 409 ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 410 ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 411 ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 412 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 413 ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 414 CHKMEMQ; 415 ierr = (*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx);CHKERRQ(ierr); 416 CHKMEMQ; 417 ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 418 ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 419 if (A != B) { 420 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 421 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 422 } 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "DMDASNESSetPicardLocal" 428 /*@C 429 DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration 430 431 Logically Collective 432 433 Input Arguments: 434 + dm - DM to associate callback with 435 . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 436 . func - local residual evaluation 437 - ctx - optional context for local residual evaluation 438 439 Calling sequence for func: 440 + info - DMDALocalInfo defining the subdomain to evaluate the residual on 441 . x - dimensional pointer to state at which to evaluate residual 442 . f - dimensional pointer to residual, write the residual here 443 - ctx - optional context passed above 444 445 Notes: The user must use 446 extern PetscErrorCode SNESPicardComputeFunction(SNES,Vec,Vec,void*); 447 extern PetscErrorCode SNESPicardComputeJacobian(SNES,Vec,Mat,Mat,MatStructure*,void*); 448 ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr); 449 in their code before calling this routine. 450 451 452 Level: beginner 453 454 .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 455 @*/ 456 PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*), 457 PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 458 { 459 PetscErrorCode ierr; 460 DMSNES sdm; 461 DMSNES_DA *dmdasnes; 462 463 PetscFunctionBegin; 464 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 465 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 466 ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 467 468 dmdasnes->residuallocalimode = imode; 469 dmdasnes->rhsplocal = func; 470 dmdasnes->jacobianplocal = jac; 471 dmdasnes->picardlocalctx = ctx; 472 473 ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 474 PetscFunctionReturn(0); 475 } 476