xref: /petsc/src/dm/dt/tests/ex11.c (revision 0faed27c4f59cfdcd00951238dac7483cf8e56f1)
1736d4f42SpierreXVI #include <petscfv.h>
2736d4f42SpierreXVI 
3736d4f42SpierreXVI static char help[] = "Test memory allocation of PetscFV arrays used in PetscFVComputeGradient";
4736d4f42SpierreXVI 
5736d4f42SpierreXVI int main(int argc, char **argv)
6736d4f42SpierreXVI {
7736d4f42SpierreXVI   PetscFV        fvm;
8736d4f42SpierreXVI   PetscInt       dim, numFaces;
9736d4f42SpierreXVI   PetscScalar    *dx, *grad;
10736d4f42SpierreXVI 
119566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, PETSC_NULL, help));
12736d4f42SpierreXVI 
13736d4f42SpierreXVI   /*
14736d4f42SpierreXVI    Working with a 2D mesh, made of triangles, and using the 2nd neighborhood
15736d4f42SpierreXVI    to reconstruct the cell gradient with a least square method, we use numFaces = 9
16736d4f42SpierreXVI    The array dx is not initialised, but it doesn't matter here
17736d4f42SpierreXVI   */
18736d4f42SpierreXVI   dim = 2;
19736d4f42SpierreXVI   numFaces = 9;
209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim * numFaces, &dx, dim * numFaces, &grad));
219566063dSJacob Faibussowitsch   PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm));
229566063dSJacob Faibussowitsch   PetscCall(PetscFVSetType(fvm, PETSCFVLEASTSQUARES));
239566063dSJacob Faibussowitsch   PetscCall(PetscFVLeastSquaresSetMaxFaces(fvm, numFaces));
24736d4f42SpierreXVI 
25736d4f42SpierreXVI   /* Issue here */
269566063dSJacob Faibussowitsch   PetscCall(PetscFVComputeGradient(fvm, numFaces, dx, grad));
27736d4f42SpierreXVI 
289566063dSJacob Faibussowitsch   PetscCall(PetscFVDestroy(&fvm));
299566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dx, grad));
309566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
314e2b1bc7SBarry Smith   return(0);
32736d4f42SpierreXVI }
33736d4f42SpierreXVI 
34736d4f42SpierreXVI /*TEST
35*0faed27cSBarry Smith 
36736d4f42SpierreXVI   test:
37736d4f42SpierreXVI     suffix: 1
38*0faed27cSBarry Smith 
39736d4f42SpierreXVI TEST*/
40