Actual source code: petscsysbase.h
1: !
2: ! Manually maintained part of the base include file for Fortran use of PETSc.
3: ! Note: This file should contain only define statements
4: !
5: #if !defined (PETSCSYSBASEDEF_H)
6: #define PETSCSYSBASEDEF_H
7: #include "petscconf.h"
8: #if defined (PETSC_HAVE_MPIUNI)
9: #include "petsc/mpiuni/mpiunifdef.h"
10: #endif
11: #include "petscversion.h"
13: !
14: #define integer8 integer(kind=C_INT64_T)
15: #define integer4 integer(kind=C_INT32_T)
16: #define integer2 integer(kind=C_INT16_T)
17: #define integer1 integer(kind=C_INT8_T)
18: #define PetscBool logical(kind=C_BOOL)
20: #if (PETSC_SIZEOF_VOID_P == 8)
21: #define PetscOffset integer8
22: #define PetscFortranAddr integer8
23: #else
24: #define PetscOffset integer4
25: #define PetscFortranAddr integer4
26: #endif
28: #if defined(PETSC_USE_64BIT_INDICES)
29: #define PetscInt integer8
30: #else
31: #define PetscInt integer4
32: #endif
33: #define PetscInt64 integer8
35: #if defined(PETSC_USE_64BIT_BLAS_INDICES)
36: #define PetscBLASInt integer8
37: #else
38: #define PetscBLASInt integer4
39: #endif
40: #define PetscCuBLASInt integer4
41: #define PetscHipBLASInt integer4
43: !
44: #define PetscSizeT integer(kind=C_SIZE_T)
45: !
46: #if defined(PETSC_USE_MPI_F08)
47: #define MPIU_Comm type(MPI_Comm)
48: #define MPIU_Group type(MPI_Group)
49: #define MPIU_Datatype type(MPI_Datatype)
50: #define MPIU_Op type(MPI_Op)
51: #define MPIU_Request type(MPI_Request)
52: #define MPIU_Status type(MPI_Status)
53: #else
54: #define MPIU_Comm integer4
55: #define MPIU_Group integer4
56: #define MPIU_Datatype integer4
57: #define MPIU_Op integer4
58: #define MPIU_Status integer4
59: #define MPIU_Request integer4
60: #endif
61: !
62: #define PetscEnum integer4
63: #define PetscVoid PetscFortranAddr
64: !
65: #define PetscFortranFloat real(kind=C_FLOAT)
66: #define PetscFortranDouble real(kind=C_DOUBLE)
67: #define PetscFortranLongDouble real(kind=C_FLOAT128)
68: #if defined(PETSC_USE_REAL_SINGLE)
69: #define PetscComplex complex(kind=C_FLOAT_COMPLEX)
70: #elif defined(PETSC_USE_REAL_DOUBLE)
71: #define PetscComplex complex(kind=C_DOUBLE_COMPLEX)
72: #elif defined(PETSC_USE_REAL___FLOAT128)
73: #define PetscComplex complex(kind=C_FLOAT128_COMPLEX)
74: #endif
76: #if defined(PETSC_USE_COMPLEX)
77: #define PETSC_SCALAR PETSC_COMPLEX
78: #else
79: #if defined(PETSC_USE_REAL_SINGLE)
80: #define PETSC_SCALAR PETSC_FLOAT
81: #elif defined(PETSC_USE_REAL___FLOAT128)
82: #define PETSC_SCALAR PETSC___FLOAT128
83: #else
84: #define PETSC_SCALAR PETSC_DOUBLE
85: #endif
86: #endif
87: #if defined(PETSC_USE_REAL_SINGLE)
88: #define PETSC_REAL PETSC_FLOAT
89: #define PetscIntToReal(a) real(a)
90: #elif defined(PETSC_USE_REAL___FLOAT128)
91: #define PETSC_REAL PETSC___FLOAT128
92: #define PetscIntToReal(a) dble(a)
93: #else
94: #define PETSC_REAL PETSC_DOUBLE
95: #define PetscIntToReal(a) dble(a)
96: #endif
97: !
98: ! Macro for templating between real and complex
99: !
100: #if defined(PETSC_USE_COMPLEX)
101: #define PetscScalar PetscComplex
102: !
103: ! F90 uses real(), conjg() when KIND parameter is used.
104: !
105: #define PetscRealPart(a) real(a)
106: #define PetscConj(a) conjg(a)
107: #define PetscImaginaryPart(a) aimag(a)
108: #else
109: #if defined (PETSC_USE_REAL_SINGLE)
110: #define PetscScalar PetscFortranFloat
111: #elif defined(PETSC_USE_REAL___FLOAT128)
112: #define PetscScalar PetscFortranLongDouble
113: #elif defined(PETSC_USE_REAL_DOUBLE)
114: #define PetscScalar PetscFortranDouble
115: #endif
116: #define PetscRealPart(a) a
117: #define PetscConj(a) a
118: #define PetscImaginaryPart(a) 0.0
119: #endif
121: #if defined (PETSC_USE_REAL_SINGLE)
122: #define PetscReal PetscFortranFloat
123: #elif defined(PETSC_USE_REAL___FLOAT128)
124: #define PetscReal PetscFortranLongDouble
125: #elif defined(PETSC_USE_REAL_DOUBLE)
126: #define PetscReal PetscFortranDouble
127: #endif
129: #define PetscReal2d type(tPetscReal2d)
131: #define PETSC_FORTRAN_TYPE_INITIALIZE -2
132: #define PetscObjectIsNull(obj) (obj%v == 0 .or. obj%v == PETSC_FORTRAN_TYPE_INITIALIZE .or. obj%v == -3)
133: #define PetscObjectNullify(obj) obj%v = PETSC_FORTRAN_TYPE_INITIALIZE
134: !
135: ! Macros for error checking
136: !
137: #define SETERRQ(c, ierr, s) call PetscError(c, ierr, PETSC_ERROR_INITIAL, s); return
138: #define SETERRA(c, ierr, s) call PetscError(c, ierr, PETSC_ERROR_INITIAL, s); call MPIU_Abort(c, ierr)
139: #if defined(PETSC_HAVE_FORTRAN_FREE_LINE_LENGTH_NONE)
140: #define CHKERRQ(ierr) if (ierr .ne. 0) then;call PetscErrorF(ierr,__LINE__,__FILE__);return;endif
141: #define CHKERRA(ierr) if (ierr .ne. 0) then;call PetscErrorF(ierr,__LINE__,__FILE__);call MPIU_Abort(PETSC_COMM_SELF,ierr);endif
142: #define CHKERRMPI(ierr) if (ierr .ne. 0) then;call PetscErrorMPI(ierr,__LINE__,__FILE__);return;endif
143: #define CHKERRMPIA(ierr) if (ierr .ne. 0) then;call PetscErrorMPI(ierr,__LINE__,__FILE__);call MPIU_Abort(PETSC_COMM_SELF,ierr);endif
144: #else
145: #define CHKERRQ(ierr) if (ierr .ne. 0) then;call PetscErrorF(ierr);return;endif
146: #define CHKERRA(ierr) if (ierr .ne. 0) then;call PetscErrorF(ierr);call MPIU_Abort(PETSC_COMM_SELF,ierr);endif
147: #define CHKERRMPI(ierr) if (ierr .ne. 0) then;call PetscErrorMPI(ierr);return;endif
148: #define CHKERRMPIA(ierr) if (ierr .ne. 0) then;call PetscErrorMPI(ierr);call MPIU_Abort(PETSC_COMM_SELF,ierr);endif
149: #endif
150: #define CHKMEMQ call chkmemfortran(__LINE__,__FILE__,ierr)
151: #define PetscCall(func) call func; CHKERRQ(ierr)
152: #define PetscCallMPI(func) call func; CHKERRMPI(ierr)
153: #define PetscCallA(func) call func; CHKERRA(ierr)
154: #define PetscCallMPIA(func) call func; CHKERRMPIA(ierr)
155: #define PetscCheckA(err, c, ierr, s) if (.not.(err)) then; SETERRA(c, ierr, s); endif
156: #define PetscCheck(err, c, ierr, s) if (.not.(err)) then; SETERRQ(c, ierr, s); endif
158: #if !defined(PetscFlush)
159: #if defined(PETSC_HAVE_FORTRAN_FLUSH)
160: #define PetscFlush(a) flush(a)
161: #elif defined(PETSC_HAVE_FORTRAN_FLUSH_)
162: #define PetscFlush(a) flush_(a)
163: #else
164: #define PetscFlush(a)
165: #endif
166: #endif
168: #define PetscEnumCase(e) case(e%v)
170: #define PetscObjectSpecificCast(sp,ob) sp%v = ob%v
172: #endif