1 #include <petsc-private/dmimpl.h> 2 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 3 4 typedef struct { 5 PetscErrorCode (*residuallocal)(DM,Vec,Vec,void*); 6 PetscErrorCode (*jacobianlocal)(DM,Vec,Mat,Mat,MatStructure*,void*); 7 void *residuallocalctx; 8 void *jacobianlocalctx; 9 } DMSNES_Local; 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "DMSNESDestroy_DMLocal" 13 static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm) 14 { 15 PetscErrorCode ierr; 16 17 PetscFunctionBegin; 18 ierr = PetscFree(sdm->data);CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 22 #undef __FUNCT__ 23 #define __FUNCT__ "DMSNESDuplicate_DMLocal" 24 static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm) 25 { 26 PetscErrorCode ierr; 27 28 PetscFunctionBegin; 29 ierr = PetscNewLog(sdm,DMSNES_Local,&sdm->data);CHKERRQ(ierr); 30 if (oldsdm->data) { 31 ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));CHKERRQ(ierr); 32 } 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "DMLocalSNESGetContext" 38 static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes) 39 { 40 PetscErrorCode ierr; 41 42 PetscFunctionBegin; 43 *dmlocalsnes = NULL; 44 if (!sdm->data) { 45 ierr = PetscNewLog(dm,DMSNES_Local,&sdm->data);CHKERRQ(ierr); 46 47 sdm->ops->destroy = DMSNESDestroy_DMLocal; 48 sdm->ops->duplicate = DMSNESDuplicate_DMLocal; 49 } 50 *dmlocalsnes = (DMSNES_Local*)sdm->data; 51 PetscFunctionReturn(0); 52 } 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "SNESComputeFunction_DMLocal" 56 static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx) 57 { 58 PetscErrorCode ierr; 59 DM dm; 60 DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 61 Vec Xloc,Floc; 62 63 PetscFunctionBegin; 64 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 65 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 66 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 67 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 68 ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 69 ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 70 ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 71 ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 72 ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 73 ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 74 CHKMEMQ; 75 ierr = (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);CHKERRQ(ierr); 76 CHKMEMQ; 77 ierr = VecZeroEntries(F);CHKERRQ(ierr); 78 ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 79 ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 80 ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 81 ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "SNESComputeJacobian_DMLocal" 87 static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 88 { 89 PetscErrorCode ierr; 90 DM dm; 91 DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 92 Vec Xloc; 93 94 PetscFunctionBegin; 95 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 96 if (dmlocalsnes->jacobianlocal) { 97 ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 98 ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 99 ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 100 ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 101 CHKMEMQ; 102 ierr = (*dmlocalsnes->jacobianlocal)(dm,Xloc,*A,*B,mstr,dmlocalsnes->jacobianlocalctx);CHKERRQ(ierr); 103 CHKMEMQ; 104 ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 105 } else { 106 MatFDColoring fdcoloring; 107 ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 108 if (!fdcoloring) { 109 ISColoring coloring; 110 111 ierr = DMCreateColoring(dm,dm->coloringtype,dm->mattype,&coloring);CHKERRQ(ierr); 112 ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr); 113 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 114 switch (dm->coloringtype) { 115 case IS_COLORING_GLOBAL: 116 ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 117 break; 118 default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 119 } 120 ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 121 ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 122 ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 123 ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 124 125 /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 126 * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 127 * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 128 * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 129 * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 130 */ 131 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 132 } 133 *mstr = SAME_NONZERO_PATTERN; 134 ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr); 135 } 136 /* This will be redundant if the user called both, but it's too common to forget. */ 137 if (*A != *B) { 138 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 139 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 140 } 141 PetscFunctionReturn(0); 142 } 143 144 #undef __FUNCT__ 145 #define __FUNCT__ "DMSNESSetFunctionLocal" 146 /*@C 147 DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 148 containing the local vector information PLUS ghost point information. It should compute a result for all local 149 elements and DMSNES will automatically accumulate the overlapping values. 150 151 Logically Collective 152 153 Input Arguments: 154 + dm - DM to associate callback with 155 . func - local residual evaluation 156 - ctx - optional context for local residual evaluation 157 158 Level: beginner 159 160 .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 161 @*/ 162 PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 163 { 164 PetscErrorCode ierr; 165 DMSNES sdm; 166 DMSNES_Local *dmlocalsnes; 167 168 PetscFunctionBegin; 169 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 170 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 171 ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 172 173 dmlocalsnes->residuallocal = func; 174 dmlocalsnes->residuallocalctx = ctx; 175 176 ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 177 if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 178 ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 179 } 180 PetscFunctionReturn(0); 181 } 182 183 #undef __FUNCT__ 184 #define __FUNCT__ "DMSNESSetJacobianLocal" 185 /*@C 186 DMSNESSetJacobianLocal - set a local Jacobian evaluation function 187 188 Logically Collective 189 190 Input Arguments: 191 + dm - DM to associate callback with 192 . func - local Jacobian evaluation 193 - ctx - optional context for local Jacobian evaluation 194 195 Level: beginner 196 197 .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 198 @*/ 199 PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,MatStructure*,void*),void *ctx) 200 { 201 PetscErrorCode ierr; 202 DMSNES sdm; 203 DMSNES_Local *dmlocalsnes; 204 205 PetscFunctionBegin; 206 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 207 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 208 ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 209 210 dmlocalsnes->jacobianlocal = func; 211 dmlocalsnes->jacobianlocalctx = ctx; 212 213 ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217