xref: /petsc/src/sys/tests/ex36f.F90 (revision f14a7c29b82d1117d8e3de344377442be395a55f)
1c4762a1bSJed Brown!
2c4762a1bSJed Brown!  Program to test PetscRandom, PetscObjectReference() and other PetscObjectXXX functions.
3c4762a1bSJed Brown!
4c4762a1bSJed Brown#include <petsc/finclude/petscsys.h>
5*c5e229c2SMartin Diehlprogram main
6c4762a1bSJed Brown  use petscsys
7c4762a1bSJed Brown  implicit none
8c4762a1bSJed Brown
9c4762a1bSJed Brown  PetscErrorCode ierr
10c4762a1bSJed Brown  PetscRandom r, q, r2
11c4762a1bSJed Brown  PetscScalar rand
12c4762a1bSJed Brown  PetscInt ref
13c4762a1bSJed Brown
14f8402805SBarry Smith  PetscCallA(PetscInitialize(ierr))
15c4762a1bSJed Brown
16f8402805SBarry Smith  PetscCallA(PetscRandomCreate(PETSC_COMM_WORLD, r, ierr))
17f8402805SBarry Smith  PetscCallA(PetscRandomCreate(PETSC_COMM_WORLD, r2, ierr))
18f8402805SBarry Smith  PetscCallA(PetscRandomSetFromOptions(r, ierr))
19f8402805SBarry Smith  PetscCallA(PetscRandomGetValue(r, rand, ierr))
20c4762a1bSJed Brown  print *, 'Random value:', rand
21c4762a1bSJed Brown
22f8402805SBarry Smith  PetscCallA(PetscObjectReference(r, ierr))
23f8402805SBarry Smith  PetscCallA(PetscObjectGetReference(r, ref, ierr))
24c4762a1bSJed Brown  print *, 'Reference value:', ref
25f8402805SBarry Smith  PetscCallA(PetscObjectDereference(r, ierr))
26c4762a1bSJed Brown
27ccfd86f1SBarry Smith  PetscCallA(PetscObjectCompose(r, 'test', r2, ierr))
28ccfd86f1SBarry Smith  PetscCallA(PetscObjectQuery(r, 'test', q, ierr))
294820e4eaSBarry Smith  PetscCheckA(q == r2, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Object compose/query failed')
30c4762a1bSJed Brown
31f8402805SBarry Smith  PetscCallA(PetscRandomDestroy(r, ierr))
32f8402805SBarry Smith  PetscCallA(PetscRandomDestroy(r2, ierr))
33f8402805SBarry Smith  PetscCallA(PetscFinalize(ierr))
34c4762a1bSJed Brownend
35c4762a1bSJed Brown
36c4762a1bSJed Brown!
37c4762a1bSJed Brown!/*TEST
38c4762a1bSJed Brown!
39c4762a1bSJed Brown!   test:
40c4762a1bSJed Brown!     requires: !complex
41c4762a1bSJed Brown!
42c4762a1bSJed Brown!TEST*/
43