Actual source code: dmlocalsnes.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/dmimpl.h>
2: #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
4: typedef struct {
5: PetscErrorCode (*residuallocal)(DM,Vec,Vec,void*);
6: PetscErrorCode (*jacobianlocal)(DM,Vec,Mat,Mat,void*);
7: void *residuallocalctx;
8: void *jacobianlocalctx;
9: } DMSNES_Local;
13: static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm)
14: {
18: PetscFree(sdm->data);
19: return(0);
20: }
24: static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm)
25: {
29: PetscNewLog(sdm,(DMSNES_Local**)&sdm->data);
30: if (oldsdm->data) {
31: PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));
32: }
33: return(0);
34: }
38: static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes)
39: {
43: *dmlocalsnes = NULL;
44: if (!sdm->data) {
45: PetscNewLog(dm,(DMSNES_Local**)&sdm->data);
47: sdm->ops->destroy = DMSNESDestroy_DMLocal;
48: sdm->ops->duplicate = DMSNESDuplicate_DMLocal;
49: }
50: *dmlocalsnes = (DMSNES_Local*)sdm->data;
51: return(0);
52: }
56: static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx)
57: {
59: DM dm;
60: DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx;
61: Vec Xloc,Floc;
67: SNESGetDM(snes,&dm);
68: DMGetLocalVector(dm,&Xloc);
69: DMGetLocalVector(dm,&Floc);
70: VecZeroEntries(Xloc);
71: VecZeroEntries(Floc);
72: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
73: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
74: CHKMEMQ;
75: (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);
76: CHKMEMQ;
77: VecZeroEntries(F);
78: DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);
79: DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);
80: DMRestoreLocalVector(dm,&Floc);
81: DMRestoreLocalVector(dm,&Xloc);
82: return(0);
83: }
87: static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat A,Mat B,void *ctx)
88: {
90: DM dm;
91: DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx;
92: Vec Xloc;
95: SNESGetDM(snes,&dm);
96: if (dmlocalsnes->jacobianlocal) {
97: DMGetLocalVector(dm,&Xloc);
98: VecZeroEntries(Xloc);
99: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
100: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
101: CHKMEMQ;
102: (*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx);
103: CHKMEMQ;
104: DMRestoreLocalVector(dm,&Xloc);
105: } else {
106: MatFDColoring fdcoloring;
107: PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);
108: if (!fdcoloring) {
109: ISColoring coloring;
111: DMCreateColoring(dm,dm->coloringtype,&coloring);
112: MatFDColoringCreate(B,coloring,&fdcoloring);
113: ISColoringDestroy(&coloring);
114: switch (dm->coloringtype) {
115: case IS_COLORING_GLOBAL:
116: MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);
117: break;
118: default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
119: }
120: PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);
121: MatFDColoringSetFromOptions(fdcoloring);
122: MatFDColoringSetUp(B,coloring,fdcoloring);
123: PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);
124: PetscObjectDereference((PetscObject)fdcoloring);
126: /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
127: * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
128: * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
129: * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
130: * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
131: */
132: PetscObjectDereference((PetscObject)dm);
133: }
134: MatFDColoringApply(B,fdcoloring,X,snes);
135: }
136: /* This will be redundant if the user called both, but it's too common to forget. */
137: if (A != B) {
138: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
139: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
140: }
141: return(0);
142: }
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.
151: Logically Collective
153: Input Arguments:
154: + dm - DM to associate callback with
155: . func - local residual evaluation
156: - ctx - optional context for local residual evaluation
158: Level: beginner
160: .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
161: @*/
162: PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx)
163: {
165: DMSNES sdm;
166: DMSNES_Local *dmlocalsnes;
170: DMGetDMSNESWrite(dm,&sdm);
171: DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);
173: dmlocalsnes->residuallocal = func;
174: dmlocalsnes->residuallocalctx = ctx;
176: DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);
177: if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
178: DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);
179: }
180: return(0);
181: }
185: /*@C
186: DMSNESSetJacobianLocal - set a local Jacobian evaluation function
188: Logically Collective
190: Input Arguments:
191: + dm - DM to associate callback with
192: . func - local Jacobian evaluation
193: - ctx - optional context for local Jacobian evaluation
195: Level: beginner
197: .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
198: @*/
199: PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,void*),void *ctx)
200: {
202: DMSNES sdm;
203: DMSNES_Local *dmlocalsnes;
207: DMGetDMSNESWrite(dm,&sdm);
208: DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);
210: dmlocalsnes->jacobianlocal = func;
211: dmlocalsnes->jacobianlocalctx = ctx;
213: DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);
214: return(0);
215: }