Actual source code: dmlocalsnes.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <petsc/private/dmimpl.h>
  2:  #include <petsc/private/snesimpl.h>

  4: typedef struct {
  5:   PetscErrorCode (*residuallocal)(DM,Vec,Vec,void*);
  6:   PetscErrorCode (*jacobianlocal)(DM,Vec,Mat,Mat,void*);
  7:   PetscErrorCode (*boundarylocal)(DM,Vec,void*);
  8:   void *residuallocalctx;
  9:   void *jacobianlocalctx;
 10:   void *boundarylocalctx;
 11: } DMSNES_Local;

 13: static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm)
 14: {

 18:   PetscFree(sdm->data);
 19:   return(0);
 20: }

 22: static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm)
 23: {

 27:   PetscNewLog(sdm,(DMSNES_Local**)&sdm->data);
 28:   if (oldsdm->data) {
 29:     PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));
 30:   }
 31:   return(0);
 32: }

 34: static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes)
 35: {

 39:   *dmlocalsnes = NULL;
 40:   if (!sdm->data) {
 41:     PetscNewLog(dm,(DMSNES_Local**)&sdm->data);

 43:     sdm->ops->destroy   = DMSNESDestroy_DMLocal;
 44:     sdm->ops->duplicate = DMSNESDuplicate_DMLocal;
 45:   }
 46:   *dmlocalsnes = (DMSNES_Local*)sdm->data;
 47:   return(0);
 48: }

 50: static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx)
 51: {
 53:   DM             dm;
 54:   DMSNES_Local   *dmlocalsnes = (DMSNES_Local*)ctx;
 55:   Vec            Xloc,Floc;

 61:   SNESGetDM(snes,&dm);
 62:   DMGetLocalVector(dm,&Xloc);
 63:   DMGetLocalVector(dm,&Floc);
 64:   VecZeroEntries(Xloc);
 65:   if (dmlocalsnes->boundarylocal) {(*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);}
 66:   VecZeroEntries(Floc);
 67:   DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
 68:   DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
 69:   CHKMEMQ;
 70:   (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);
 71:   CHKMEMQ;
 72:   VecZeroEntries(F);
 73:   DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);
 74:   DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);
 75:   DMRestoreLocalVector(dm,&Floc);
 76:   DMRestoreLocalVector(dm,&Xloc);
 77:   {
 78:     char        name[PETSC_MAX_PATH_LEN];
 79:     char        oldname[PETSC_MAX_PATH_LEN];
 80:     const char *tmp;
 81:     PetscInt    it;

 83:     SNESGetIterationNumber(snes, &it);
 84:     PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int) it);
 85:     PetscObjectGetName((PetscObject) X, &tmp);
 86:     PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN-1);
 87:     PetscObjectSetName((PetscObject) X, name);
 88:     VecViewFromOptions(X, (PetscObject) snes, "-dmsnes_solution_vec_view");
 89:     PetscObjectSetName((PetscObject) X, oldname);
 90:     PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int) it);
 91:     PetscObjectSetName((PetscObject) F, name);
 92:     VecViewFromOptions(F, (PetscObject) snes, "-dmsnes_residual_vec_view");
 93:   }
 94:   return(0);
 95: }

 97: static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat A,Mat B,void *ctx)
 98: {
100:   DM             dm;
101:   DMSNES_Local   *dmlocalsnes = (DMSNES_Local*)ctx;
102:   Vec            Xloc;

105:   SNESGetDM(snes,&dm);
106:   if (dmlocalsnes->jacobianlocal) {
107:     DMGetLocalVector(dm,&Xloc);
108:     VecZeroEntries(Xloc);
109:     if (dmlocalsnes->boundarylocal) {(*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);}
110:     DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
111:     DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
112:     CHKMEMQ;
113:     (*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx);
114:     CHKMEMQ;
115:     DMRestoreLocalVector(dm,&Xloc);
116:   } else {
117:     MatFDColoring fdcoloring;
118:     PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);
119:     if (!fdcoloring) {
120:       ISColoring coloring;

122:       DMCreateColoring(dm,dm->coloringtype,&coloring);
123:       MatFDColoringCreate(B,coloring,&fdcoloring);
124:       ISColoringDestroy(&coloring);
125:       switch (dm->coloringtype) {
126:       case IS_COLORING_GLOBAL:
127:         MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);
128:         break;
129:       default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
130:       }
131:       PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);
132:       MatFDColoringSetFromOptions(fdcoloring);
133:       MatFDColoringSetUp(B,coloring,fdcoloring);
134:       PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);
135:       PetscObjectDereference((PetscObject)fdcoloring);

137:       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
138:        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
139:        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
140:        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
141:        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
142:        */
143:       PetscObjectDereference((PetscObject)dm);
144:     }
145:     MatFDColoringApply(B,fdcoloring,X,snes);
146:   }
147:   /* This will be redundant if the user called both, but it's too common to forget. */
148:   if (A != B) {
149:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
150:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
151:   }
152:   return(0);
153: }

155: /*@C
156:    DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector
157:       containing the local vector information PLUS ghost point information. It should compute a result for all local
158:       elements and DMSNES will automatically accumulate the overlapping values.

160:    Logically Collective

162:    Input Arguments:
163: +  dm - DM to associate callback with
164: .  func - local residual evaluation
165: -  ctx - optional context for local residual evaluation

167:    Level: beginner

169: .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
170: @*/
171: PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx)
172: {
174:   DMSNES         sdm;
175:   DMSNES_Local   *dmlocalsnes;

179:   DMGetDMSNESWrite(dm,&sdm);
180:   DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);

182:   dmlocalsnes->residuallocal    = func;
183:   dmlocalsnes->residuallocalctx = ctx;

185:   DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);
186:   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
187:     DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);
188:   }
189:   return(0);
190: }

192: /*@C
193:    DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector
194:       containing the local vector information PLUS ghost point information. It should insert values into the local
195:       vector that do not come from the global vector, such as essential boundary condition data.

197:    Logically Collective

199:    Input Arguments:
200: +  dm - DM to associate callback with
201: .  func - local boundary value evaluation
202: -  ctx - optional context for local boundary value evaluation

204:    Level: intermediate

206: .seealso: DMSNESSetFunctionLocal(), DMDASNESSetJacobianLocal()
207: @*/
208: PetscErrorCode DMSNESSetBoundaryLocal(DM dm,PetscErrorCode (*func)(DM,Vec,void*),void *ctx)
209: {
211:   DMSNES         sdm;
212:   DMSNES_Local   *dmlocalsnes;

216:   DMGetDMSNESWrite(dm,&sdm);
217:   DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);

219:   dmlocalsnes->boundarylocal    = func;
220:   dmlocalsnes->boundarylocalctx = ctx;

222:   return(0);
223: }

225: /*@C
226:    DMSNESSetJacobianLocal - set a local Jacobian evaluation function

228:    Logically Collective

230:    Input Arguments:
231: +  dm - DM to associate callback with
232: .  func - local Jacobian evaluation
233: -  ctx - optional context for local Jacobian evaluation

235:    Level: beginner

237: .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
238: @*/
239: PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,void*),void *ctx)
240: {
242:   DMSNES         sdm;
243:   DMSNES_Local   *dmlocalsnes;

247:   DMGetDMSNESWrite(dm,&sdm);
248:   DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);

250:   dmlocalsnes->jacobianlocal    = func;
251:   dmlocalsnes->jacobianlocalctx = ctx;

253:   DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);
254:   return(0);
255: }