1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdm.h> /*I "petscdm.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "PetscContainerDestroy_SNESDM" 6 static PetscErrorCode PetscContainerDestroy_SNESDM(void *ctx) 7 { 8 PetscErrorCode ierr; 9 SNESDM sdm = (SNESDM)ctx; 10 11 PetscFunctionBegin; 12 if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);} 13 ierr = PetscFree(sdm);CHKERRQ(ierr); 14 PetscFunctionReturn(0); 15 } 16 17 #undef __FUNCT__ 18 #define __FUNCT__ "DMCoarsenHook_SNESDM" 19 /* Attaches the SNESDM to the coarse level. 20 * Under what conditions should we copy versus duplicate? 21 */ 22 static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx) 23 { 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr); 28 PetscFunctionReturn(0); 29 } 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "DMRestrictHook_SNESDM" 33 /* This could restrict auxiliary information to the coarse level. 34 */ 35 static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 36 { 37 38 PetscFunctionBegin; 39 PetscFunctionReturn(0); 40 } 41 42 #undef __FUNCT__ 43 #define __FUNCT__ "DMSNESGetContext" 44 /*@C 45 DMSNESGetContext - get read-only private SNESDM context from a DM 46 47 Not Collective 48 49 Input Argument: 50 . dm - DM to be used with SNES 51 52 Output Argument: 53 . snesdm - private SNESDM context 54 55 Level: developer 56 57 Notes: 58 Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible. 59 60 .seealso: DMSNESGetContextWrite() 61 @*/ 62 PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm) 63 { 64 PetscErrorCode ierr; 65 PetscContainer container; 66 SNESDM sdm; 67 68 PetscFunctionBegin; 69 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 70 ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr); 71 if (container) { 72 ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr); 73 } else { 74 ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr); 75 ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 76 ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr); 77 ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 78 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr); 79 ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 80 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 81 ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr); 82 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 83 } 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "DMSNESGetContextWrite" 89 /*@C 90 DMSNESGetContextWrite - get write access to private SNESDM context from a DM 91 92 Not Collective 93 94 Input Argument: 95 . dm - DM to be used with SNES 96 97 Output Argument: 98 . snesdm - private SNESDM context 99 100 Level: developer 101 102 .seealso: DMSNESGetContext() 103 @*/ 104 PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm) 105 { 106 PetscErrorCode ierr; 107 SNESDM sdm; 108 109 PetscFunctionBegin; 110 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 111 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 112 if (!sdm->originaldm) sdm->originaldm = dm; 113 if (sdm->originaldm != dm) { /* Copy on write */ 114 PetscContainer container; 115 SNESDM oldsdm = sdm; 116 ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr); 117 ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 118 ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr); 119 ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr); 120 ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 121 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr); 122 ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 123 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 124 } 125 *snesdm = sdm; 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "DMSNESCopyContext" 131 /*@C 132 DMSNESCopyContext - copies a DM context to a new DM 133 134 Logically Collective 135 136 Input Arguments: 137 + dmsrc - DM to obtain context from 138 - dmdest - DM to add context to 139 140 Level: developer 141 142 Note: 143 The context is copied by reference. This function does not ensure that a context exists. 144 145 .seealso: DMSNESGetContext(), SNESSetDM() 146 @*/ 147 PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest) 148 { 149 PetscErrorCode ierr; 150 PetscContainer container; 151 152 PetscFunctionBegin; 153 PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 154 PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 155 ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr); 156 if (container) { 157 ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 158 ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 159 } 160 PetscFunctionReturn(0); 161 } 162 163 #undef __FUNCT__ 164 #define __FUNCT__ "DMSNESSetFunction" 165 /*@C 166 DMSNESSetFunction - set SNES residual evaluation function 167 168 Not Collective 169 170 Input Arguments: 171 + dm - DM to be used with SNES 172 . func - residual evaluation function, see SNESSetFunction() for calling sequence 173 - ctx - context for residual evaluation 174 175 Level: advanced 176 177 Note: 178 SNESSetFunction() is normally used, but it calls this function internally because the user context is actually 179 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 180 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 181 182 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 183 @*/ 184 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 185 { 186 PetscErrorCode ierr; 187 SNESDM sdm; 188 189 PetscFunctionBegin; 190 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 191 ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 192 if (func) sdm->computefunction = func; 193 if (ctx) sdm->functionctx = ctx; 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "DMSNESGetFunction" 199 /*@C 200 DMSNESGetFunction - get SNES residual evaluation function 201 202 Not Collective 203 204 Input Argument: 205 . dm - DM to be used with SNES 206 207 Output Arguments: 208 + func - residual evaluation function, see SNESSetFunction() for calling sequence 209 - ctx - context for residual evaluation 210 211 Level: advanced 212 213 Note: 214 SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 215 associated with the DM. 216 217 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction() 218 @*/ 219 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 220 { 221 PetscErrorCode ierr; 222 SNESDM sdm; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 226 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 227 if (func) *func = sdm->computefunction; 228 if (ctx) *ctx = sdm->functionctx; 229 PetscFunctionReturn(0); 230 } 231 232 #undef __FUNCT__ 233 #define __FUNCT__ "DMSNESSetGS" 234 /*@C 235 DMSNESSetGS - set SNES Gauss-Seidel relaxation function 236 237 Not Collective 238 239 Input Argument: 240 + dm - DM to be used with SNES 241 . func - relaxation function, see SNESSetGS() for calling sequence 242 - ctx - context for residual evaluation 243 244 Level: advanced 245 246 Note: 247 SNESSetGS() is normally used, but it calls this function internally because the user context is actually 248 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 249 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 250 251 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction() 252 @*/ 253 PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 254 { 255 PetscErrorCode ierr; 256 SNESDM sdm; 257 258 PetscFunctionBegin; 259 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 260 ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 261 if (func) sdm->computegs = func; 262 if (ctx) sdm->gsctx = ctx; 263 PetscFunctionReturn(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "DMSNESGetGS" 268 /*@C 269 DMSNESGetGS - get SNES Gauss-Seidel relaxation function 270 271 Not Collective 272 273 Input Argument: 274 . dm - DM to be used with SNES 275 276 Output Arguments: 277 + func - relaxation function, see SNESSetGS() for calling sequence 278 - ctx - context for residual evaluation 279 280 Level: advanced 281 282 Note: 283 SNESGetGS() is normally used, but it calls this function internally because the user context is actually 284 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 285 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 286 287 .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction() 288 @*/ 289 PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 290 { 291 PetscErrorCode ierr; 292 SNESDM sdm; 293 294 PetscFunctionBegin; 295 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 296 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 297 if (func) *func = sdm->computegs; 298 if (ctx) *ctx = sdm->gsctx; 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "DMSNESSetJacobian" 304 /*@C 305 DMSNESSetFunction - set SNES Jacobian evaluation function 306 307 Not Collective 308 309 Input Argument: 310 + dm - DM to be used with SNES 311 . func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 312 - ctx - context for residual evaluation 313 314 Level: advanced 315 316 Note: 317 SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually 318 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 319 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 320 321 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian() 322 @*/ 323 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 324 { 325 PetscErrorCode ierr; 326 SNESDM sdm; 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 330 ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 331 if (func) sdm->computejacobian = func; 332 if (ctx) sdm->jacobianctx = ctx; 333 PetscFunctionReturn(0); 334 } 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "DMSNESGetJacobian" 338 /*@C 339 DMSNESGetFunction - get SNES Jacobian evaluation function 340 341 Not Collective 342 343 Input Argument: 344 . dm - DM to be used with SNES 345 346 Output Arguments: 347 + func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 348 - ctx - context for residual evaluation 349 350 Level: advanced 351 352 Note: 353 SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually 354 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 355 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 356 357 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 358 @*/ 359 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 360 { 361 PetscErrorCode ierr; 362 SNESDM sdm; 363 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 366 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 367 if (func) *func = sdm->computejacobian; 368 if (ctx) *ctx = sdm->jacobianctx; 369 PetscFunctionReturn(0); 370 } 371 372 #undef __FUNCT__ 373 #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy" 374 static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx) 375 { 376 PetscErrorCode ierr; 377 DM dm; 378 379 PetscFunctionBegin; 380 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 381 ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy" 387 static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 388 { 389 PetscErrorCode ierr; 390 DM dm; 391 392 PetscFunctionBegin; 393 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 394 ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "DMSNESSetUpLegacy" 400 /* Sets up calling of legacy DM routines */ 401 PetscErrorCode DMSNESSetUpLegacy(DM dm) 402 { 403 PetscErrorCode ierr; 404 SNESDM sdm; 405 406 PetscFunctionBegin; 407 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 408 if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 409 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 410 if (!sdm->computejacobian) {ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 411 PetscFunctionReturn(0); 412 } 413