xref: /petsc/src/snes/utils/dmdasnes.c (revision 6cab3a1b8d4fca83e1d7b61474c608460de73d5f)
1*6cab3a1bSJed Brown #include <petscdmda.h>          /*I "petscdmda.h" I*/
2*6cab3a1bSJed Brown #include <private/snesimpl.h>   /*I "petscsnes.h" I*/
3*6cab3a1bSJed Brown 
4*6cab3a1bSJed Brown /* This structure holds the user-provided DMDA callbacks */
5*6cab3a1bSJed Brown typedef struct {
6*6cab3a1bSJed Brown   PetscErrorCode (*residuallocal)(DMDALocalInfo*,void*,void*,void*);
7*6cab3a1bSJed Brown   PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*);
8*6cab3a1bSJed Brown   void *residuallocalctx;
9*6cab3a1bSJed Brown   void *jacobianlocalctx;
10*6cab3a1bSJed Brown } DM_DA_SNES;
11*6cab3a1bSJed Brown 
12*6cab3a1bSJed Brown #undef __FUNCT__
13*6cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA"
14*6cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM sdm)
15*6cab3a1bSJed Brown {
16*6cab3a1bSJed Brown   PetscErrorCode ierr;
17*6cab3a1bSJed Brown 
18*6cab3a1bSJed Brown   PetscFunctionBegin;
19*6cab3a1bSJed Brown   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
20*6cab3a1bSJed Brown   PetscFunctionReturn(0);
21*6cab3a1bSJed Brown }
22*6cab3a1bSJed Brown 
23*6cab3a1bSJed Brown #undef __FUNCT__
24*6cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext"
25*6cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes)
26*6cab3a1bSJed Brown {
27*6cab3a1bSJed Brown   PetscErrorCode ierr;
28*6cab3a1bSJed Brown 
29*6cab3a1bSJed Brown   PetscFunctionBegin;
30*6cab3a1bSJed Brown   if (!sdm->data) {
31*6cab3a1bSJed Brown     ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr);
32*6cab3a1bSJed Brown     sdm->destroy = SNESDMDestroy_DMDA;
33*6cab3a1bSJed Brown   }
34*6cab3a1bSJed Brown   *dmdasnes = (DM_DA_SNES*)sdm->data;
35*6cab3a1bSJed Brown   PetscFunctionReturn(0);
36*6cab3a1bSJed Brown }
37*6cab3a1bSJed Brown 
38*6cab3a1bSJed Brown #undef __FUNCT__
39*6cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA"
40*6cab3a1bSJed Brown /*
41*6cab3a1bSJed Brown   This function should eventually replace:
42*6cab3a1bSJed Brown     DMDAComputeFunction() and DMDAComputeFunction1()
43*6cab3a1bSJed Brown  */
44*6cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
45*6cab3a1bSJed Brown {
46*6cab3a1bSJed Brown   PetscErrorCode ierr;
47*6cab3a1bSJed Brown   DM             dm;
48*6cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
49*6cab3a1bSJed Brown   DMDALocalInfo  info;
50*6cab3a1bSJed Brown   Vec            Xloc;
51*6cab3a1bSJed Brown   void           *x,*f;
52*6cab3a1bSJed Brown 
53*6cab3a1bSJed Brown   PetscFunctionBegin;
54*6cab3a1bSJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
55*6cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
56*6cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
57*6cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
58*6cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
59*6cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
60*6cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
61*6cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
62*6cab3a1bSJed Brown   CHKMEMQ;
63*6cab3a1bSJed Brown   ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
64*6cab3a1bSJed Brown   CHKMEMQ;
65*6cab3a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
66*6cab3a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
67*6cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
68*6cab3a1bSJed Brown   PetscFunctionReturn(0);
69*6cab3a1bSJed Brown }
70*6cab3a1bSJed Brown 
71*6cab3a1bSJed Brown #undef __FUNCT__
72*6cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal"
73*6cab3a1bSJed Brown /*@C
74*6cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
75*6cab3a1bSJed Brown 
76*6cab3a1bSJed Brown    Logically Collective
77*6cab3a1bSJed Brown 
78*6cab3a1bSJed Brown    Input Arguments:
79*6cab3a1bSJed Brown +  dm - DM to associate callback with
80*6cab3a1bSJed Brown .  func - local residual evaluation
81*6cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
82*6cab3a1bSJed Brown 
83*6cab3a1bSJed Brown    Calling sequence for func:
84*6cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
85*6cab3a1bSJed Brown .  x - dimensional pointer to state at which to evaluate residual
86*6cab3a1bSJed Brown .  f - dimensional pointer to residual, write the residual here
87*6cab3a1bSJed Brown -  ctx - optional context passed above
88*6cab3a1bSJed Brown 
89*6cab3a1bSJed Brown    Level: beginner
90*6cab3a1bSJed Brown 
91*6cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
92*6cab3a1bSJed Brown @*/
93*6cab3a1bSJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
94*6cab3a1bSJed Brown {
95*6cab3a1bSJed Brown   PetscErrorCode ierr;
96*6cab3a1bSJed Brown   SNESDM         sdm;
97*6cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes;
98*6cab3a1bSJed Brown 
99*6cab3a1bSJed Brown   PetscFunctionBegin;
100*6cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
101*6cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
102*6cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
103*6cab3a1bSJed Brown   dmdasnes->residuallocal = func;
104*6cab3a1bSJed Brown   dmdasnes->residuallocalctx = ctx;
105*6cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
106*6cab3a1bSJed Brown   PetscFunctionReturn(0);
107*6cab3a1bSJed Brown }
108