xref: /petsc/src/sys/tests/ex10.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests PetscArraymove()/PetscMemmove()\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscsys.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   PetscInt       i,*a,*b;
9c4762a1bSJed Brown 
10*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
11c4762a1bSJed Brown 
125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(10,&a));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(20,&b));
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   /*
16c4762a1bSJed Brown       Nonoverlapping regions
17c4762a1bSJed Brown   */
18c4762a1bSJed Brown   for (i=0; i<20; i++) b[i] = i;
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraymove(a,b,10));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(10,a,NULL));
21c4762a1bSJed Brown 
225f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraymove(a,b,15));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(15,a,NULL));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(b));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /*
35c4762a1bSJed Brown      |       |                    |       |
36c4762a1bSJed Brown      a       b                   a+20   a+25
37c4762a1bSJed Brown                                         b+20
38c4762a1bSJed Brown   */
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(25,&a));
40c4762a1bSJed Brown   b    = a + 5;
41c4762a1bSJed Brown   for (i=0; i<20; i++) b[i] = i;
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraymove(a,b,20));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscIntView(20,a,NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(a));
45c4762a1bSJed Brown 
46*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
47*b122ec5aSJacob Faibussowitsch   return 0;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown /*TEST
51c4762a1bSJed Brown 
52c4762a1bSJed Brown    test:
53c4762a1bSJed Brown 
54c4762a1bSJed Brown TEST*/
55