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