xref: /petsc/src/ts/utils/dmplexts.c (revision b59f628eff0cd0541d0457de20841eec6c6bc428)
1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2 #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
3 #include <petscds.h>
4 #include <petscfv.h>
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "DMPlexTSGetGeometryFVM"
8 /*@
9   DMPlexTSGetGeometryFVM - Return precomputed geometric data
10 
11   Input Parameter:
12 . dm - The DM
13 
14   Output Parameters:
15 + facegeom - The values precomputed from face geometry
16 . cellgeom - The values precomputed from cell geometry
17 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
18 
19   Level: developer
20 
21 .seealso: DMPlexTSSetRHSFunctionLocal()
22 @*/
23 PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
24 {
25   DMTS           dmts;
26   PetscObject    obj;
27   PetscErrorCode ierr;
28 
29   PetscFunctionBegin;
30   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
31   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
32   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
33   if (!obj) {
34     Vec cellgeom, facegeom;
35 
36     ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr);
37     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
38     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
39     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
40     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
41   }
42   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
43   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
44   if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);}
45   PetscFunctionReturn(0);
46 }
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "DMPlexTSGetGradientDM"
50 /*@C
51   DMPlexTSGetGradientDM - Return gradient data layout
52 
53   Input Parameters:
54 + dm - The DM
55 - fv - The PetscFV
56 
57   Output Parameter:
58 . dmGrad - The layout for gradient values
59 
60   Level: developer
61 
62 .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
63 @*/
64 PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
65 {
66   DMTS           dmts;
67   PetscObject    obj;
68   PetscBool      computeGradients;
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
73   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
74   PetscValidPointer(dmGrad,3);
75   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
76   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
77   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
78   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
79   if (!obj) {
80     DM  dmGrad;
81     Vec faceGeometry, cellGeometry;
82 
83     ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
84     ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
85     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
86     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
87   }
88   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
94 /*@
95   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
96 
97   Input Parameters:
98 + dm - The mesh
99 . t - The time
100 . locX  - Local solution
101 - user - The user context
102 
103   Output Parameter:
104 . locF  - Local output vector
105 
106   Level: developer
107 
108 .seealso: DMPlexComputeJacobianActionFEM()
109 @*/
110 PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec locF, void *user)
111 {
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   ierr = DMPlexComputeResidual_Internal(dm, time, locX, NULL, locF, user);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
121 /*@
122   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
123 
124   Input Parameters:
125 + dm - The mesh
126 . t - The time
127 . locX  - Local solution
128 . locX_t - Local solution time derivative, or NULL
129 - user - The user context
130 
131   Output Parameter:
132 . locF  - Local output vector
133 
134   Level: developer
135 
136 .seealso: DMPlexComputeJacobianActionFEM()
137 @*/
138 PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
139 {
140   PetscErrorCode ierr;
141 
142   PetscFunctionBegin;
143   ierr = DMPlexComputeResidual_Internal(dm, time, locX, locX_t, locF, user);CHKERRQ(ierr);
144   PetscFunctionReturn(0);
145 }
146