xref: /petsc/src/sys/tests/ex10.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Tests PetscArraymove()/PetscMemmove()\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsys.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   PetscInt i, *a, *b;
8c4762a1bSJed Brown 
9327415f7SBarry Smith   PetscFunctionBeginUser;
10*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
11c4762a1bSJed Brown 
129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(10, &a));
139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(20, &b));
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   /*
16c4762a1bSJed Brown       Nonoverlapping regions
17c4762a1bSJed Brown   */
18c4762a1bSJed Brown   for (i = 0; i < 20; i++) b[i] = i;
199566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(a, b, 10));
209566063dSJacob Faibussowitsch   PetscCall(PetscIntView(10, a, NULL));
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(PetscFree(a));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /*
25c4762a1bSJed Brown      |        |                |       |
26c4762a1bSJed Brown      b        a               b+15    b+20
27c4762a1bSJed Brown                               a+10    a+15
28c4762a1bSJed Brown   */
29c4762a1bSJed Brown   a = b + 5;
309566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(a, b, 15));
319566063dSJacob Faibussowitsch   PetscCall(PetscIntView(15, a, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /*
35c4762a1bSJed Brown      |       |                    |       |
36c4762a1bSJed Brown      a       b                   a+20   a+25
37c4762a1bSJed Brown                                         b+20
38c4762a1bSJed Brown   */
399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(25, &a));
40c4762a1bSJed Brown   b = a + 5;
41c4762a1bSJed Brown   for (i = 0; i < 20; i++) b[i] = i;
429566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(a, b, 20));
439566063dSJacob Faibussowitsch   PetscCall(PetscIntView(20, a, NULL));
449566063dSJacob Faibussowitsch   PetscCall(PetscFree(a));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
47b122ec5aSJacob Faibussowitsch   return 0;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown /*TEST
51c4762a1bSJed Brown 
52c4762a1bSJed Brown    test:
53c4762a1bSJed Brown 
54c4762a1bSJed Brown TEST*/
55