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