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