1736d4f42SpierreXVI #include <petscfv.h>
2736d4f42SpierreXVI
3736d4f42SpierreXVI static char help[] = "Test memory allocation of PetscFV arrays used in PetscFVComputeGradient";
4736d4f42SpierreXVI
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7736d4f42SpierreXVI PetscFV fvm;
8736d4f42SpierreXVI PetscInt dim, numFaces;
9736d4f42SpierreXVI PetscScalar *dx, *grad;
10736d4f42SpierreXVI
11327415f7SBarry Smith PetscFunctionBeginUser;
12f3fa974cSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
13736d4f42SpierreXVI
14736d4f42SpierreXVI /*
15736d4f42SpierreXVI Working with a 2D mesh, made of triangles, and using the 2nd neighborhood
16736d4f42SpierreXVI to reconstruct the cell gradient with a least square method, we use numFaces = 9
17736d4f42SpierreXVI The array dx is not initialised, but it doesn't matter here
18736d4f42SpierreXVI */
19736d4f42SpierreXVI dim = 2;
20736d4f42SpierreXVI numFaces = 9;
219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dim * numFaces, &dx, dim * numFaces, &grad));
229566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm));
239566063dSJacob Faibussowitsch PetscCall(PetscFVSetType(fvm, PETSCFVLEASTSQUARES));
249566063dSJacob Faibussowitsch PetscCall(PetscFVLeastSquaresSetMaxFaces(fvm, numFaces));
25736d4f42SpierreXVI
26736d4f42SpierreXVI /* Issue here */
279566063dSJacob Faibussowitsch PetscCall(PetscFVComputeGradient(fvm, numFaces, dx, grad));
28736d4f42SpierreXVI
299566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fvm));
309566063dSJacob Faibussowitsch PetscCall(PetscFree2(dx, grad));
319566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
324ad8454bSPierre Jolivet return 0;
33736d4f42SpierreXVI }
34736d4f42SpierreXVI
35736d4f42SpierreXVI /*TEST
360faed27cSBarry Smith
37736d4f42SpierreXVI test:
38736d4f42SpierreXVI suffix: 1
39*3886731fSPierre Jolivet output_file: output/empty.out
400faed27cSBarry Smith
41736d4f42SpierreXVI TEST*/
42