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