xref: /petsc/src/ts/utils/dmplexts.c (revision 3b16df4212dd91db9c2bec74fb927e3c2384b1bd)
11faf85eaSMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2254c1ad2SMatthew G. Knepley #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
3924a1b8fSMatthew G. Knepley #include <petscds.h>
46dbbd306SMatthew G. Knepley #include <petscfv.h>
56dbbd306SMatthew G. Knepley 
6254c1ad2SMatthew G. Knepley #undef __FUNCT__
7b6aca0f9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometryFVM"
8c510411aSMatthew G. Knepley /*@
9b6aca0f9SMatthew G. Knepley   DMPlexTSGetGeometryFVM - Return precomputed geometric data
10c510411aSMatthew G. Knepley 
11c510411aSMatthew G. Knepley   Input Parameter:
12c510411aSMatthew G. Knepley . dm - The DM
13c510411aSMatthew G. Knepley 
14c510411aSMatthew G. Knepley   Output Parameters:
15c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
16c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
17c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
18c510411aSMatthew G. Knepley 
19c510411aSMatthew G. Knepley   Level: developer
20c510411aSMatthew G. Knepley 
21c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
22c510411aSMatthew G. Knepley @*/
23b6aca0f9SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
24a0ac79e7SMatthew G. Knepley {
25a0ac79e7SMatthew G. Knepley   DMTS           dmts;
261faf85eaSMatthew G. Knepley   PetscObject    obj;
27a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
28a0ac79e7SMatthew G. Knepley 
29a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
30924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
31a0ac79e7SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
32b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
331faf85eaSMatthew G. Knepley   if (!obj) {
341faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
351faf85eaSMatthew G. Knepley 
361faf85eaSMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr);
37b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
38b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
391faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
401faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
411faf85eaSMatthew G. Knepley   }
42b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
43b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
44113c68e6SMatthew G. Knepley   if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);}
45924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
46924a1b8fSMatthew G. Knepley }
47924a1b8fSMatthew G. Knepley 
48924a1b8fSMatthew G. Knepley #undef __FUNCT__
49c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
50c49ccbb3SMatthew G. Knepley /*@C
51c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
52c49ccbb3SMatthew G. Knepley 
53c49ccbb3SMatthew G. Knepley   Input Parameters:
54c49ccbb3SMatthew G. Knepley + dm - The DM
55c49ccbb3SMatthew G. Knepley - fv - The PetscFV
56c49ccbb3SMatthew G. Knepley 
57c49ccbb3SMatthew G. Knepley   Output Parameter:
58c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
59c49ccbb3SMatthew G. Knepley 
60c49ccbb3SMatthew G. Knepley   Level: developer
61c49ccbb3SMatthew G. Knepley 
62c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
63c49ccbb3SMatthew G. Knepley @*/
64c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
65c49ccbb3SMatthew G. Knepley {
66c49ccbb3SMatthew G. Knepley   DMTS           dmts;
67c49ccbb3SMatthew G. Knepley   PetscObject    obj;
68c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
69c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
70c49ccbb3SMatthew G. Knepley 
71c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
72c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
73c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
74c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
75c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
76c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
77c49ccbb3SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
78c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
79c49ccbb3SMatthew G. Knepley   if (!obj) {
80c49ccbb3SMatthew G. Knepley     DM  dmGrad;
8108449791SMatthew G. Knepley     Vec faceGeometry, cellGeometry;
82c49ccbb3SMatthew G. Knepley 
8308449791SMatthew G. Knepley     ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
8408449791SMatthew G. Knepley     ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
85c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
86924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
87c49ccbb3SMatthew G. Knepley   }
88c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
89254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
90254c1ad2SMatthew G. Knepley }
91254c1ad2SMatthew G. Knepley 
92254c1ad2SMatthew G. Knepley #undef __FUNCT__
93924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
9408449791SMatthew G. Knepley /*@
9508449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
9608449791SMatthew G. Knepley 
9708449791SMatthew G. Knepley   Input Parameters:
9808449791SMatthew G. Knepley + dm - The mesh
9908449791SMatthew G. Knepley . t - The time
10008449791SMatthew G. Knepley . locX  - Local solution
10108449791SMatthew G. Knepley - user - The user context
10208449791SMatthew G. Knepley 
10308449791SMatthew G. Knepley   Output Parameter:
104*3b16df42SMatthew G. Knepley . F  - Global output vector
10508449791SMatthew G. Knepley 
10608449791SMatthew G. Knepley   Level: developer
10708449791SMatthew G. Knepley 
10808449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
10908449791SMatthew G. Knepley @*/
110*3b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
111254c1ad2SMatthew G. Knepley {
112*3b16df42SMatthew G. Knepley   Vec            locF;
113254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
114254c1ad2SMatthew G. Knepley 
115254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
116*3b16df42SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locF);CHKERRQ(ierr);
117*3b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
11808449791SMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(dm, time, locX, NULL, locF, user);CHKERRQ(ierr);
119*3b16df42SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, locF, INSERT_VALUES, F);CHKERRQ(ierr);
120*3b16df42SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, locF, INSERT_VALUES, F);CHKERRQ(ierr);
121*3b16df42SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locF);CHKERRQ(ierr);
122254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
123254c1ad2SMatthew G. Knepley }
124254c1ad2SMatthew G. Knepley 
125254c1ad2SMatthew G. Knepley #undef __FUNCT__
12624cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
12724cdb843SMatthew G. Knepley /*@
12824cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
12924cdb843SMatthew G. Knepley 
13024cdb843SMatthew G. Knepley   Input Parameters:
13124cdb843SMatthew G. Knepley + dm - The mesh
13224cdb843SMatthew G. Knepley . t - The time
13308449791SMatthew G. Knepley . locX  - Local solution
13408449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
13524cdb843SMatthew G. Knepley - user - The user context
13624cdb843SMatthew G. Knepley 
13724cdb843SMatthew G. Knepley   Output Parameter:
13808449791SMatthew G. Knepley . locF  - Local output vector
13924cdb843SMatthew G. Knepley 
14024cdb843SMatthew G. Knepley   Level: developer
14124cdb843SMatthew G. Knepley 
14224cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
14324cdb843SMatthew G. Knepley @*/
14408449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
14524cdb843SMatthew G. Knepley {
14624cdb843SMatthew G. Knepley   PetscErrorCode ierr;
14724cdb843SMatthew G. Knepley 
14824cdb843SMatthew G. Knepley   PetscFunctionBegin;
14908449791SMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(dm, time, locX, locX_t, locF, user);CHKERRQ(ierr);
15024cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
15124cdb843SMatthew G. Knepley }
1527cdb2a12SMatthew G. Knepley 
1537cdb2a12SMatthew G. Knepley #undef __FUNCT__
1547cdb2a12SMatthew G. Knepley #define __FUNCT__ "DMTSCheckFromOptions"
1557871a739SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u, void (**exactFuncs)(const PetscReal x[], PetscScalar *u, void *ctx), void *ctx)
1567cdb2a12SMatthew G. Knepley {
1577cdb2a12SMatthew G. Knepley   DM             dm;
1587cdb2a12SMatthew G. Knepley   SNES           snes;
1597cdb2a12SMatthew G. Knepley   Mat            J, M;
1607cdb2a12SMatthew G. Knepley   Vec            sol, r, b;
1617cdb2a12SMatthew G. Knepley   MatNullSpace   nullSpace;
162cd45c955SMatthew G. Knepley   PetscReal     *error, res = 0.0;
163cd45c955SMatthew G. Knepley   PetscInt       numFields;
1647cdb2a12SMatthew G. Knepley   PetscBool      check;
1657cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
1667cdb2a12SMatthew G. Knepley 
1677cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
1687cdb2a12SMatthew G. Knepley   ierr = PetscOptionsHasName(ts->hdr.prefix, "-dmts_check", &check);CHKERRQ(ierr);
1697cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
1707cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
1717cdb2a12SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
1727cdb2a12SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1737cdb2a12SMatthew G. Knepley   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
1747cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
1757cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
1767cdb2a12SMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
1777cdb2a12SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
1787cdb2a12SMatthew G. Knepley   M    = J;
1797cdb2a12SMatthew G. Knepley #if 0
1807cdb2a12SMatthew G. Knepley   {
1817cdb2a12SMatthew G. Knepley     ierr = CreatePressureNullSpace(dm, user, &nullSpace);CHKERRQ(ierr);
1827cdb2a12SMatthew G. Knepley     ierr = MatSetNullSpace(J, nullSpace);CHKERRQ(ierr);
1837cdb2a12SMatthew G. Knepley     if (M != J) {ierr = MatSetNullSpace(M, nullSpace);CHKERRQ(ierr);}
1847cdb2a12SMatthew G. Knepley   }
1857cdb2a12SMatthew G. Knepley #endif
1867cdb2a12SMatthew G. Knepley   /* Check discretization error */
187cd45c955SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
188cd45c955SMatthew G. Knepley   ierr = PetscMalloc1(PetscMax(1, numFields), &error);CHKERRQ(ierr);
189cd45c955SMatthew G. Knepley   if (numFields > 1) {
190cd45c955SMatthew G. Knepley     PetscInt f;
191cd45c955SMatthew G. Knepley 
1927871a739SMatthew G. Knepley     ierr = DMPlexComputeL2FieldDiff(dm, exactFuncs, ctx, u, error);CHKERRQ(ierr);
193cd45c955SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: [");CHKERRQ(ierr);
194cd45c955SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
195cd45c955SMatthew G. Knepley       if (f) {ierr = PetscPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);}
196cd45c955SMatthew G. Knepley       if (error[f] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "%g", error[f]);CHKERRQ(ierr);}
197cd45c955SMatthew G. Knepley       else                     {ierr = PetscPrintf(PETSC_COMM_WORLD, "< 1.0e-11");CHKERRQ(ierr);}
198cd45c955SMatthew G. Knepley     }
199cd45c955SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "]\n");CHKERRQ(ierr);
200cd45c955SMatthew G. Knepley   } else {
2017871a739SMatthew G. Knepley     ierr = DMPlexComputeL2Diff(dm, exactFuncs, ctx, u, &error[0]);CHKERRQ(ierr);
202cd45c955SMatthew G. Knepley     if (error[0] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error[0]);CHKERRQ(ierr);}
203cd45c955SMatthew G. Knepley     else                     {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
204cd45c955SMatthew G. Knepley   }
205cd45c955SMatthew G. Knepley   ierr = PetscFree(error);CHKERRQ(ierr);
2067cdb2a12SMatthew G. Knepley   /* Check residual */
2077cdb2a12SMatthew G. Knepley   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
2087cdb2a12SMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2097cdb2a12SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);CHKERRQ(ierr);
2107cdb2a12SMatthew G. Knepley   ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2117cdb2a12SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
2127cdb2a12SMatthew G. Knepley   ierr = VecViewFromOptions(r, "res_", "-vec_view");CHKERRQ(ierr);
2137cdb2a12SMatthew G. Knepley   /* Check Jacobian */
2147cdb2a12SMatthew G. Knepley   ierr = SNESComputeJacobian(snes, u, M, M);CHKERRQ(ierr);
2157cdb2a12SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullSpace);CHKERRQ(ierr);
2167cdb2a12SMatthew G. Knepley   if (nullSpace) {
2177cdb2a12SMatthew G. Knepley     PetscBool isNull;
2187cdb2a12SMatthew G. Knepley     ierr = MatNullSpaceTest(nullSpace, J, &isNull);CHKERRQ(ierr);
2197cdb2a12SMatthew G. Knepley     if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
2207cdb2a12SMatthew G. Knepley   }
2217cdb2a12SMatthew G. Knepley   ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
2227cdb2a12SMatthew G. Knepley   ierr = VecSet(r, 0.0);CHKERRQ(ierr);
2237cdb2a12SMatthew G. Knepley   ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
2247cdb2a12SMatthew G. Knepley   ierr = MatMult(M, u, r);CHKERRQ(ierr);
2257cdb2a12SMatthew G. Knepley   ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
2267cdb2a12SMatthew G. Knepley   ierr = VecDestroy(&b);CHKERRQ(ierr);
2277cdb2a12SMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2287cdb2a12SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);CHKERRQ(ierr);
2297cdb2a12SMatthew G. Knepley   ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2307cdb2a12SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) r, "Au - b = Au + F(0)");CHKERRQ(ierr);
2317cdb2a12SMatthew G. Knepley   ierr = VecViewFromOptions(r, "linear_res_", "-vec_view");CHKERRQ(ierr);
2327cdb2a12SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2337cdb2a12SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
2347cdb2a12SMatthew G. Knepley   ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
2357cdb2a12SMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
2367cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
2377cdb2a12SMatthew G. Knepley }
238