Actual source code: f90_fwrap.F90

  1: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
  3: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  4: #include <petsc/finclude/petscsys.h>
  5: subroutine F90Array1dCreateScalar(array, start, len1, ptr)
  6:   use, intrinsic :: ISO_C_binding
  7:   implicit none
  8:   PetscInt start, len1
  9:   PetscScalar, target :: array(start:start + len1 - 1)
 10:   PetscScalar, pointer :: ptr(:)

 12:   ptr => array
 13: end subroutine

 15: subroutine F90Array1dCreateReal(array, start, len1, ptr)
 16:   use, intrinsic :: ISO_C_binding
 17:   implicit none
 18:   PetscInt start, len1
 19:   PetscReal, target :: array(start:start + len1 - 1)
 20:   PetscReal, pointer :: ptr(:)

 22:   ptr => array
 23: end subroutine

 25: subroutine F90Array1dCreateInt(array, start, len1, ptr)
 26:   use, intrinsic :: ISO_C_binding
 27:   implicit none
 28:   PetscInt start, len1
 29:   PetscInt, target :: array(start:start + len1 - 1)
 30:   PetscInt, pointer :: ptr(:)

 32:   ptr => array
 33: end subroutine

 35: subroutine F90Array1dCreateMPIInt(array, start, len1, ptr)
 36:   use, intrinsic :: ISO_C_binding
 37:   implicit none
 38:   PetscInt start, len1
 39:   PetscMPIInt, target :: array(start:start + len1 - 1)
 40:   PetscMPIInt, pointer :: ptr(:)

 42:   ptr => array
 43: end subroutine

 45: subroutine F90Array1dCreateFortranAddr(array, start, len1, ptr)
 46:   use, intrinsic :: ISO_C_binding
 47:   implicit none
 48:   PetscInt start, len1
 49:   PetscFortranAddr, target :: array(start:start + len1 - 1)
 50:   PetscFortranAddr, pointer :: ptr(:)

 52:   ptr => array
 53: end subroutine

 55: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 56: subroutine F90Array1dAccessScalar(ptr, address)
 57:   use, intrinsic :: ISO_C_binding
 58:   implicit none
 59:   PetscScalar, pointer :: ptr(:)
 60:   PetscFortranAddr address
 61:   PetscInt start

 63:   if (.not. associated(ptr)) then
 64:     address = 0
 65:   else
 66:     start = lbound(ptr, 1)
 67:     call F90Array1dGetAddrScalar(ptr(start), address)
 68:   end if
 69: end subroutine

 71: subroutine F90Array1dAccessReal(ptr, address)
 72:   use, intrinsic :: ISO_C_binding
 73:   implicit none
 74:   PetscReal, pointer :: ptr(:)
 75:   PetscFortranAddr address
 76:   PetscInt start

 78:   if (.not. associated(ptr)) then
 79:     address = 0
 80:   else
 81:     start = lbound(ptr, 1)
 82:     call F90Array1dGetAddrReal(ptr(start), address)
 83:   end if
 84: end subroutine

 86: subroutine F90Array1dAccessInt(ptr, address)
 87:   use, intrinsic :: ISO_C_binding
 88:   implicit none
 89:   PetscInt, pointer :: ptr(:)
 90:   PetscFortranAddr address
 91:   PetscInt start

 93:   if (.not. associated(ptr)) then
 94:     address = 0
 95:   else
 96:     start = lbound(ptr, 1)
 97:     call F90Array1dGetAddrInt(ptr(start), address)
 98:   end if
 99: end subroutine

101: subroutine F90Array1dAccessMPIInt(ptr, address)
102:   use, intrinsic :: ISO_C_binding
103:   implicit none
104:   PetscMPIInt, pointer :: ptr(:)
105:   PetscFortranAddr address
106:   PetscInt start

108:   if (.not. associated(ptr)) then
109:     address = 0
110:   else
111:     start = lbound(ptr, 1)
112:     call F90Array1dGetAddrMPIInt(ptr(start), address)
113:   end if
114: end subroutine

116: subroutine F90Array1dAccessFortranAddr(ptr, address)
117:   use, intrinsic :: ISO_C_binding
118:   implicit none
119:   PetscFortranAddr, pointer :: ptr(:)
120:   PetscFortranAddr address
121:   PetscInt start

123:   if (.not. associated(ptr)) then
124:     address = 0
125:   else
126:     start = lbound(ptr, 1)
127:     call F90Array1dGetAddrFortranAddr(ptr(start), address)
128:   end if
129: end subroutine

131: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
132: subroutine F90Array1dDestroyScalar(ptr)
133:   use, intrinsic :: ISO_C_binding
134:   implicit none
135:   PetscScalar, pointer :: ptr(:)

137:   nullify (ptr)
138: end subroutine

140: subroutine F90Array1dDestroyReal(ptr)
141:   use, intrinsic :: ISO_C_binding
142:   implicit none
143:   PetscReal, pointer :: ptr(:)

145:   nullify (ptr)
146: end subroutine

148: subroutine F90Array1dDestroyInt(ptr)
149:   use, intrinsic :: ISO_C_binding
150:   implicit none
151:   PetscInt, pointer :: ptr(:)

153:   nullify (ptr)
154: end subroutine

156: subroutine F90Array1dDestroyMPIInt(ptr)
157:   use, intrinsic :: ISO_C_binding
158:   implicit none
159:   PetscMPIInt, pointer :: ptr(:)

161:   nullify (ptr)
162: end subroutine

164: subroutine F90Array1dDestroyFortranAddr(ptr)
165:   use, intrinsic :: ISO_C_binding
166:   implicit none
167:   PetscFortranAddr, pointer :: ptr(:)

169:   nullify (ptr)
170: end subroutine
171: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
172: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
173: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
174: subroutine F90Array2dCreateScalar(array, start1, len1, start2, len2, ptr)
175:   use, intrinsic :: ISO_C_binding
176:   implicit none
177:   PetscInt start1, len1
178:   PetscInt start2, len2
179:   PetscScalar, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1)
180:   PetscScalar, pointer :: ptr(:, :)

182:   ptr => array
183: end subroutine

185: subroutine F90Array2dCreateReal(array, start1, len1, start2, len2, ptr)
186:   use, intrinsic :: ISO_C_binding
187:   implicit none
188:   PetscInt start1, len1
189:   PetscInt start2, len2
190:   PetscReal, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1)
191:   PetscReal, pointer :: ptr(:, :)

193:   ptr => array
194: end subroutine

196: subroutine F90Array2dCreateInt(array, start1, len1, start2, len2, ptr)
197:   use, intrinsic :: ISO_C_binding
198:   implicit none
199:   PetscInt start1, len1
200:   PetscInt start2, len2
201:   PetscInt, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1)
202:   PetscInt, pointer :: ptr(:, :)

204:   ptr => array
205: end subroutine

207: subroutine F90Array2dCreateFortranAddr(array, start1, len1, start2, len2, ptr)
208:   use, intrinsic :: ISO_C_binding
209:   implicit none
210:   PetscInt start1, len1
211:   PetscInt start2, len2
212:   PetscFortranAddr, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1)
213:   PetscFortranAddr, pointer :: ptr(:, :)

215:   ptr => array
216: end subroutine

218: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
219: subroutine F90Array2dAccessScalar(ptr, address)
220:   use, intrinsic :: ISO_C_binding
221:   implicit none
222:   PetscScalar, pointer :: ptr(:, :)
223:   PetscFortranAddr address
224:   PetscInt start1, start2

226:   start1 = lbound(ptr, 1)
227:   start2 = lbound(ptr, 2)
228:   call F90Array2dGetAddrScalar(ptr(start1, start2), address)
229: end subroutine

231: subroutine F90Array2dAccessReal(ptr, address)
232:   use, intrinsic :: ISO_C_binding
233:   implicit none
234:   PetscReal, pointer :: ptr(:, :)
235:   PetscFortranAddr address
236:   PetscInt start1, start2

238:   start1 = lbound(ptr, 1)
239:   start2 = lbound(ptr, 2)
240:   call F90Array2dGetAddrReal(ptr(start1, start2), address)
241: end subroutine

243: subroutine F90Array2dAccessInt(ptr, address)
244:   use, intrinsic :: ISO_C_binding
245:   implicit none
246:   PetscInt, pointer :: ptr(:, :)
247:   PetscFortranAddr address
248:   PetscInt start1, start2

250:   start1 = lbound(ptr, 1)
251:   start2 = lbound(ptr, 2)
252:   call F90Array2dGetAddrInt(ptr(start1, start2), address)
253: end subroutine

255: subroutine F90Array2dAccessFortranAddr(ptr, address)
256:   use, intrinsic :: ISO_C_binding
257:   implicit none
258:   PetscFortranAddr, pointer :: ptr(:, :)
259:   PetscFortranAddr address
260:   PetscInt start1, start2

262:   start1 = lbound(ptr, 1)
263:   start2 = lbound(ptr, 2)
264:   call F90Array2dGetAddrFortranAddr(ptr(start1, start2), address)
265: end subroutine

267: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
268: subroutine F90Array2dDestroyScalar(ptr)
269:   use, intrinsic :: ISO_C_binding
270:   implicit none
271:   PetscScalar, pointer :: ptr(:, :)

273:   nullify (ptr)
274: end subroutine

276: subroutine F90Array2dDestroyReal(ptr)
277:   use, intrinsic :: ISO_C_binding
278:   implicit none
279:   PetscReal, pointer :: ptr(:, :)

281:   nullify (ptr)
282: end subroutine

284: subroutine F90Array2dDestroyInt(ptr)
285:   use, intrinsic :: ISO_C_binding
286:   implicit none
287:   PetscInt, pointer :: ptr(:, :)

289:   nullify (ptr)
290: end subroutine

292: subroutine F90Array2dDestroyFortranAddr(ptr)
293:   use, intrinsic :: ISO_C_binding
294:   implicit none
295:   PetscFortranAddr, pointer :: ptr(:, :)

297:   nullify (ptr)
298: end subroutine
299: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
300: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
301: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
302: subroutine F90Array3dCreateScalar(array, start1, len1, start2, len2, start3, len3, ptr)
303:   use, intrinsic :: ISO_C_binding
304:   implicit none
305:   PetscInt start1, len1
306:   PetscInt start2, len2
307:   PetscInt start3, len3
308:   PetscScalar, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1)
309:   PetscScalar, pointer :: ptr(:, :, :)

311:   ptr => array
312: end subroutine

314: subroutine F90Array3dCreateReal(array, start1, len1, start2, len2, start3, len3, ptr)
315:   use, intrinsic :: ISO_C_binding
316:   implicit none
317:   PetscInt start1, len1
318:   PetscInt start2, len2
319:   PetscInt start3, len3
320:   PetscReal, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1)
321:   PetscReal, pointer :: ptr(:, :, :)

323:   ptr => array
324: end subroutine

326: subroutine F90Array3dCreateInt(array, start1, len1, start2, len2, start3, len3, ptr)
327:   use, intrinsic :: ISO_C_binding
328:   implicit none
329:   PetscInt start1, len1
330:   PetscInt start2, len2
331:   PetscInt start3, len3
332:   PetscInt, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1)
333:   PetscInt, pointer :: ptr(:, :, :)

335:   ptr => array
336: end subroutine

338: subroutine F90Array3dCreateFortranAddr(array, start1, len1, start2, len2, start3, len3, ptr)
339:   use, intrinsic :: ISO_C_binding
340:   implicit none
341:   PetscInt start1, len1
342:   PetscInt start2, len2
343:   PetscInt start3, len3
344:   PetscFortranAddr, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1)
345:   PetscFortranAddr, pointer :: ptr(:, :, :)

347:   ptr => array
348: end subroutine

350: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
351: subroutine F90Array3dAccessScalar(ptr, address)
352:   use, intrinsic :: ISO_C_binding
353:   implicit none
354:   PetscScalar, pointer :: ptr(:, :, :)
355:   PetscFortranAddr address
356:   PetscInt start1, start2, start3

358:   start1 = lbound(ptr, 1)
359:   start2 = lbound(ptr, 2)
360:   start3 = lbound(ptr, 3)
361:   call F90Array3dGetAddrScalar(ptr(start1, start2, start3), address)
362: end subroutine

364: subroutine F90Array3dAccessReal(ptr, address)
365:   use, intrinsic :: ISO_C_binding
366:   implicit none
367:   PetscReal, pointer :: ptr(:, :, :)
368:   PetscFortranAddr address
369:   PetscInt start1, start2, start3

371:   start1 = lbound(ptr, 1)
372:   start2 = lbound(ptr, 2)
373:   start3 = lbound(ptr, 3)
374:   call F90Array3dGetAddrReal(ptr(start1, start2, start3), address)
375: end subroutine

377: subroutine F90Array3dAccessInt(ptr, address)
378:   use, intrinsic :: ISO_C_binding
379:   implicit none
380:   PetscInt, pointer :: ptr(:, :, :)
381:   PetscFortranAddr address
382:   PetscInt start1, start2, start3

384:   start1 = lbound(ptr, 1)
385:   start2 = lbound(ptr, 2)
386:   start3 = lbound(ptr, 3)
387:   call F90Array3dGetAddrInt(ptr(start1, start2, start3), address)
388: end subroutine

390: subroutine F90Array3dAccessFortranAddr(ptr, address)
391:   use, intrinsic :: ISO_C_binding
392:   implicit none
393:   PetscFortranAddr, pointer :: ptr(:, :, :)
394:   PetscFortranAddr address
395:   PetscInt start1, start2, start3

397:   start1 = lbound(ptr, 1)
398:   start2 = lbound(ptr, 2)
399:   start3 = lbound(ptr, 3)
400:   call F90Array3dGetAddrFortranAddr(ptr(start1, start2, start3), address)
401: end subroutine

403: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
404: subroutine F90Array3dDestroyScalar(ptr)
405:   use, intrinsic :: ISO_C_binding
406:   implicit none
407:   PetscScalar, pointer :: ptr(:, :, :)

409:   nullify (ptr)
410: end subroutine

412: subroutine F90Array3dDestroyReal(ptr)
413:   use, intrinsic :: ISO_C_binding
414:   implicit none
415:   PetscReal, pointer :: ptr(:, :, :)

417:   nullify (ptr)
418: end subroutine

420: subroutine F90Array3dDestroyInt(ptr)
421:   use, intrinsic :: ISO_C_binding
422:   implicit none
423:   PetscInt, pointer :: ptr(:, :, :)

425:   nullify (ptr)
426: end subroutine

428: subroutine F90Array3dDestroyFortranAddr(ptr)
429:   use, intrinsic :: ISO_C_binding
430:   implicit none
431:   PetscFortranAddr, pointer :: ptr(:, :, :)

433:   nullify (ptr)
434: end subroutine

436: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
437: subroutine F90Array4dCreateScalar(array, start1, len1, start2, len2, start3, len3, start4, len4, ptr)
438:   use, intrinsic :: ISO_C_binding
439:   implicit none
440:   PetscInt start1, len1
441:   PetscInt start2, len2
442:   PetscInt start3, len3
443:   PetscInt start4, len4
444:   PetscScalar, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1, start4:start4 + len4 - 1)
445:   PetscScalar, pointer :: ptr(:, :, :, :)

447:   ptr => array
448: end subroutine

450: subroutine F90Array4dCreateReal(array, start1, len1, start2, len2, start3, len3, start4, len4, ptr)
451:   use, intrinsic :: ISO_C_binding
452:   implicit none
453:   PetscInt start1, len1
454:   PetscInt start2, len2
455:   PetscInt start3, len3
456:   PetscInt start4, len4
457:   PetscReal, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1, start4:start4 + len4 - 1)
458:   PetscReal, pointer :: ptr(:, :, :, :)

460:   ptr => array
461: end subroutine

463: subroutine F90Array4dCreateInt(array, start1, len1, start2, len2, start3, len3, start4, len4, ptr)
464:   use, intrinsic :: ISO_C_binding
465:   implicit none
466:   PetscInt start1, len1
467:   PetscInt start2, len2
468:   PetscInt start3, len3
469:   PetscInt start4, len4
470:   PetscInt, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1, start4:start4 + len4 - 1)
471:   PetscInt, pointer :: ptr(:, :, :, :)

473:   ptr => array
474: end subroutine

476: subroutine F90Array4dCreateFortranAddr(array, start1, len1, start2, len2, start3, len3, start4, len4, ptr)
477:   use, intrinsic :: ISO_C_binding
478:   implicit none
479:   PetscInt start1, len1
480:   PetscInt start2, len2
481:   PetscInt start3, len3
482:   PetscInt start4, len4
483:   PetscFortranAddr, target :: array(start1:start1 + len1 - 1, start2:start2 + len2 - 1, start3:start3 + len3 - 1, start4:start4 + len4 - 1)
484:   PetscFortranAddr, pointer :: ptr(:, :, :, :)

486:   ptr => array
487: end subroutine

489: subroutine F90Array4dAccessScalar(ptr, address)
490:   use, intrinsic :: ISO_C_binding
491:   implicit none
492:   PetscScalar, pointer :: ptr(:, :, :, :)
493:   PetscFortranAddr address
494:   PetscInt start1, start2, start3, start4

496:   start1 = lbound(ptr, 1)
497:   start2 = lbound(ptr, 2)
498:   start3 = lbound(ptr, 3)
499:   start4 = lbound(ptr, 4)
500:   call F90Array4dGetAddrScalar(ptr(start1, start2, start3, start4), address)
501: end subroutine

503: subroutine F90Array4dAccessReal(ptr, address)
504:   use, intrinsic :: ISO_C_binding
505:   implicit none
506:   PetscReal, pointer :: ptr(:, :, :, :)
507:   PetscFortranAddr address
508:   PetscInt start1, start2, start3, start4

510:   start1 = lbound(ptr, 1)
511:   start2 = lbound(ptr, 2)
512:   start3 = lbound(ptr, 3)
513:   start4 = lbound(ptr, 4)
514:   call F90Array4dGetAddrReal(ptr(start1, start2, start3, start4), address)
515: end subroutine

517: subroutine F90Array4dAccessInt(ptr, address)
518:   use, intrinsic :: ISO_C_binding
519:   implicit none
520:   PetscInt, pointer :: ptr(:, :, :, :)
521:   PetscFortranAddr address
522:   PetscInt start1, start2, start3, start4

524:   start1 = lbound(ptr, 1)
525:   start2 = lbound(ptr, 2)
526:   start3 = lbound(ptr, 3)
527:   start4 = lbound(ptr, 4)
528:   call F90Array4dGetAddrInt(ptr(start1, start2, start3, start4), address)
529: end subroutine

531: subroutine F90Array4dAccessFortranAddr(ptr, address)
532:   use, intrinsic :: ISO_C_binding
533:   implicit none
534:   PetscScalar, pointer :: ptr(:, :, :, :)
535:   PetscFortranAddr address
536:   PetscFortranAddr start1, start2, start3, start4

538:   start1 = lbound(ptr, 1)
539:   start2 = lbound(ptr, 2)
540:   start3 = lbound(ptr, 3)
541:   start4 = lbound(ptr, 4)
542:   call F90Array4dGetAddrFortranAddr(ptr(start1, start2, start3, start4), address)
543: end subroutine

545: subroutine F90Array4dDestroyScalar(ptr)
546:   use, intrinsic :: ISO_C_binding
547:   implicit none
548:   PetscScalar, pointer :: ptr(:, :, :, :)

550:   nullify (ptr)
551: end subroutine

553: subroutine F90Array4dDestroyReal(ptr)
554:   use, intrinsic :: ISO_C_binding
555:   implicit none
556:   PetscReal, pointer :: ptr(:, :, :, :)

558:   nullify (ptr)
559: end subroutine

561: subroutine F90Array4dDestroyInt(ptr)
562:   use, intrinsic :: ISO_C_binding
563:   implicit none
564:   PetscInt, pointer :: ptr(:, :, :, :)

566:   nullify (ptr)
567: end subroutine

569: subroutine F90Array4dDestroyFortranAddr(ptr)
570:   use, intrinsic :: ISO_C_binding
571:   implicit none
572:   PetscFortranAddr, pointer :: ptr(:, :, :, :)

574:   nullify (ptr)
575: end subroutine

577: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
578: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
579: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!