Actual source code: ex6f90.F
petsc-dev 2014-02-02
1: !-----------------------------------------------------------------------
2: !
3: ! manages tokamak edge region.
4: !
5: ! This is a translation of ex6.c into F90 by Alan Glasser, LANL
6: !-----------------------------------------------------------------------
7: !-----------------------------------------------------------------------
8: ! code organization.
9: !-----------------------------------------------------------------------
10: ! 0. barry_mod.
11: ! 1. barry_get_global_corners.
12: ! 2. barry_get_global_vector.
13: ! 3. barry_get_local_vector.
14: ! 4. barry_global_to_local.
15: ! 5. barry_destroy_fa.
16: ! 6. barry_create_fa.
17: ! 7. barry_draw_patch.
18: ! 8. barry_draw_fa.
19: ! 9. barry_map_region3.
20: ! 10. barry_map_region2.
21: ! 11. barry_map_region1.
22: ! 12. barry_main.
23: !-----------------------------------------------------------------------
24: ! subprogram 0. barry_mod.
25: ! module declarations.
26: !-----------------------------------------------------------------------
27: !-----------------------------------------------------------------------
28: ! declarations.
29: !-----------------------------------------------------------------------
30: MODULE barry_mod
31: IMPLICIT NONE
33: #include <finclude/petscsys.h>
34: #include <finclude/petscvec.h>
35: #include <finclude/petscdmda.h>
36: #include <finclude/petscmat.h>
37: #include <finclude/petscis.h>
38: #include <finclude/petscao.h>
39: #include <finclude/petscviewer.h>
40: #include <finclude/petscdraw.h>
42: TYPE :: fa_type
43: MPI_Comm, DIMENSION(0:2) :: comm
44: INTEGER, DIMENSION(0:2) :: xl,yl,ml,nl,xg,yg,mg,ng,offl,offg
45: Vec :: g,l
46: VecScatter :: vscat
47: INTEGER :: p1,p2,r1,r2,r1g,r2g,sw
48: END TYPE fa_type
50: TYPE :: patch_type
51: INTEGER :: mx,my
52: PetscScalar, DIMENSION(:,:,:), POINTER :: xy
53: END TYPE patch_type
55: LOGICAL, PARAMETER :: diagnose=.FALSE.
56: INTEGER :: ierr
57: PetscReal,PARAMETER::pi=3.1415926535897932385_8
58: PetscReal,Parameter::twopi=2*pi
60: CONTAINS
61: !-----------------------------------------------------------------------
62: ! subprogram 1. barry_get_global_corners.
63: ! returns global corner data.
64: !-----------------------------------------------------------------------
65: !-----------------------------------------------------------------------
66: ! declarations.
67: !-----------------------------------------------------------------------
68: SUBROUTINE barry_get_global_corners(fa,j,x,y,m,n)
70: TYPE(fa_type), INTENT(IN) :: fa
71: INTEGER, INTENT(IN) :: j
72: INTEGER, INTENT(OUT) :: x,y,m,n
73: !-----------------------------------------------------------------------
74: ! computations.
75: !-----------------------------------------------------------------------
76: IF (fa%comm(j) /= 0) THEN
77: x=fa%xg(j)
78: y=fa%yg(j)
79: m=fa%mg(j)
80: n=fa%ng(j)
81: ELSE
82: x=0
83: y=0
84: m=0
85: n=0
86: ENDIF
87: !-----------------------------------------------------------------------
88: ! terminate.
89: !-----------------------------------------------------------------------
90: RETURN
91: END SUBROUTINE barry_get_global_corners
92: !-----------------------------------------------------------------------
93: ! subprogram 2. barry_get_global_vector.
94: ! returns local vector.
95: !-----------------------------------------------------------------------
96: !-----------------------------------------------------------------------
97: ! declarations.
98: !-----------------------------------------------------------------------
99: SUBROUTINE barry_get_global_vector(fa,v)
101: TYPE(fa_type), INTENT(IN) :: fa
102: Vec, INTENT(OUT) :: v
103: !-----------------------------------------------------------------------
104: ! computations.
105: !-----------------------------------------------------------------------
106: CALL VecDuplicate(fa%g,v,ierr)
107: !-----------------------------------------------------------------------
108: ! terminate.
109: !-----------------------------------------------------------------------
110: RETURN
111: END SUBROUTINE barry_get_global_vector
112: !-----------------------------------------------------------------------
113: ! subprogram 3. barry_get_local_vector.
114: ! returns local vector.
115: !-----------------------------------------------------------------------
116: !-----------------------------------------------------------------------
117: ! declarations.
118: !-----------------------------------------------------------------------
119: SUBROUTINE barry_get_local_vector(fa,v)
121: TYPE(fa_type), INTENT(IN) :: fa
122: Vec, INTENT(OUT) :: v
123: !-----------------------------------------------------------------------
124: ! computations.
125: !-----------------------------------------------------------------------
126: CALL VecDuplicate(fa%l,v,ierr)
127: !-----------------------------------------------------------------------
128: ! terminate.
129: !-----------------------------------------------------------------------
130: RETURN
131: END SUBROUTINE barry_get_local_vector
132: !-----------------------------------------------------------------------
133: ! subprogram 4. barry_global_to_local.
134: ! performs VecScatter.
135: !-----------------------------------------------------------------------
136: !-----------------------------------------------------------------------
137: ! declarations.
138: !-----------------------------------------------------------------------
139: SUBROUTINE barry_global_to_local(fa,g,l)
141: TYPE(fa_type), INTENT(IN) :: fa
142: Vec, INTENT(IN) :: g
143: Vec, INTENT(OUT) :: l
144: !-----------------------------------------------------------------------
145: ! computations.
146: !-----------------------------------------------------------------------
147: CALL VecScatterBegin(fa%vscat,g,l,INSERT_VALUES,SCATTER_FORWARD, &
148: & ierr)
149: CALL VecScatterEnd(fa%vscat,g,l,INSERT_VALUES,SCATTER_FORWARD, &
150: & ierr)
151: !-----------------------------------------------------------------------
152: ! terminate.
153: !-----------------------------------------------------------------------
154: RETURN
155: END SUBROUTINE barry_global_to_local
156: !-----------------------------------------------------------------------
157: ! subprogram 5. barry_destroy_fa.
158: ! destroys fa.
159: !-----------------------------------------------------------------------
160: !-----------------------------------------------------------------------
161: ! declarations.
162: !-----------------------------------------------------------------------
163: SUBROUTINE barry_destroy_fa(fa)
165: TYPE(fa_type), INTENT(OUT) :: fa
166: !-----------------------------------------------------------------------
167: ! computations.
168: !-----------------------------------------------------------------------
169: CALL VecDestroy(fa%g,ierr)
170: CALL VecDestroy(fa%l,ierr)
171: CALL VecScatterDestroy(fa%vscat,ierr)
172: !-----------------------------------------------------------------------
173: ! terminate.
174: !-----------------------------------------------------------------------
175: RETURN
176: END SUBROUTINE barry_destroy_fa
177: !-----------------------------------------------------------------------
178: ! subprogram 6. barry_create_fa.
179: ! creates fa.
180: !-----------------------------------------------------------------------
181: !-----------------------------------------------------------------------
182: ! declarations.
183: !-----------------------------------------------------------------------
184: SUBROUTINE barry_create_fa(fa)
186: TYPE(fa_type), INTENT(OUT) :: fa
188: INTEGER :: rank,tonglobal,globalrstart,x,nx,y,ny,i,j,k,ig, &
189: & fromnglobal,nscat,nlocal,cntl1,cntl2,cntl3,il,it
190: INTEGER, DIMENSION(0:2) :: offt
191: INTEGER, DIMENSION(:), POINTER :: tonatural,fromnatural, &
192: & to,from,indices
193: PetscScalar, DIMENSION(1) :: globalarray
194: PetscScalar, DIMENSION(1) :: localarray
195: PetscScalar, DIMENSION(1) :: toarray
197: DM :: da1=0,da2=0,da3=0
198: Vec :: vl1=0,vl2=0,vl3=0
199: Vec :: vg1=0,vg2=0,vg3=0
200: AO :: toao,globalao
201: IS :: tois,globalis,is
202: Vec :: tovec,globalvec,localvec
203: VecScatter :: vscat
204: !-----------------------------------------------------------------------
205: ! set integer members of fa.
206: !-----------------------------------------------------------------------
207: fa%p1=24
208: fa%p2=15
209: fa%r1=6
210: fa%r2=6
211: fa%sw=1
212: fa%r1g=fa%r1+fa%sw
213: fa%r2g=fa%r2+fa%sw
214: !-----------------------------------------------------------------------
215: ! set up communicators.
216: !-----------------------------------------------------------------------
217: CALL MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
218: fa%comm=PETSC_COMM_WORLD
219: !-----------------------------------------------------------------------
220: ! set up distributed arrays.
221: !-----------------------------------------------------------------------
222: IF (fa%comm(1) /= 0) THEN
223: CALL DMDACreate2d(fa%comm(1),DMDA_BOUNDARY_PERIODIC, &
224: & DMDA_BOUNDARY_NONE, DMDA_STENCIL_BOX, &
225: & fa%p2,fa%r2g,PETSC_DECIDE,PETSC_DECIDE,1,fa%sw, &
226: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da2,ierr)
227: CALL DMGetLocalVector(da2,vl2,ierr)
228: CALL DMGetGlobalVector(da2,vg2,ierr)
229: ENDIF
230: IF (fa%comm(2) /= 0) THEN
231: CALL DMDACreate2d(fa%comm(2),DMDA_BOUNDARY_NONE, &
232: & DMDA_BOUNDARY_NONE, DMDA_STENCIL_BOX, &
233: & fa%p1-fa%p2,fa%r2g,PETSC_DECIDE,PETSC_DECIDE,1,fa%sw, &
234: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da3,ierr)
235: CALL DMGetLocalVector(da3,vl3,ierr)
236: CALL DMGetGlobalVector(da3,vg3,ierr)
237: ENDIF
238: IF (fa%comm(0) /= 0) THEN
239: CALL DMDACreate2d(fa%comm(0),DMDA_BOUNDARY_NONE, &
240: & DMDA_BOUNDARY_NONE, DMDA_STENCIL_BOX, &
241: & fa%p1,fa%r1g,PETSC_DECIDE,PETSC_DECIDE,1,fa%sw, &
242: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da1,ierr)
243: CALL DMGetLocalVector(da1,vl1,ierr)
244: CALL DMGetGlobalVector(da1,vg1,ierr)
245: ENDIF
246: !-----------------------------------------------------------------------
247: ! count the number of unknowns owned on each processor and determine
248: ! the starting point of each processors ownership
249: ! for global vector with redundancy.
250: !-----------------------------------------------------------------------
251: tonglobal = 0
252: IF (fa%comm(1) /= 0) THEN
253: CALL DMDAGetCorners(da2,x,y,0,nx,ny,0,ierr)
254: tonglobal=tonglobal+nx*ny
255: ENDIF
256: IF (fa%comm(2) /= 0) THEN
257: CALL DMDAGetCorners(da3,x,y,0,nx,ny,0,ierr)
258: tonglobal=tonglobal+nx*ny
259: ENDIF
260: IF (fa%comm(0) /= 0) THEN
261: CALL DMDAGetCorners(da1,x,y,0,nx,ny,0,ierr)
262: tonglobal=tonglobal+nx*ny
263: ENDIF
264: WRITE(*,'("[",i1,"]",a,i3)') &
265: & rank," Number of unknowns owned ",tonglobal
266: !-----------------------------------------------------------------------
267: ! get tonatural number for each node
268: !-----------------------------------------------------------------------
269: ALLOCATE(tonatural(0:tonglobal))
270: tonglobal=0
271: IF (fa%comm(1) /= 0) THEN
272: CALL DMDAGetCorners(da2,x,y,0,nx,ny,0,ierr)
273: DO j=0,ny-1
274: DO i=0,nx-1
275: tonatural(tonglobal)=(fa%p1-fa%p2)/2+x+i+fa%p1*(y+j)
276: tonglobal=tonglobal+1
277: ENDDO
278: ENDDO
279: ENDIF
280: IF (fa%comm(2) /= 0) THEN
281: CALL DMDAGetCorners(da3,x,y,0,nx,ny,0,ierr)
282: DO j=0,ny-1
283: DO i=0,nx-1
284: IF (x+i < (fa%p1-fa%p2)/2) THEN
285: tonatural(tonglobal)=x+i+fa%p1*(y+j)
286: ELSE
287: tonatural(tonglobal)=fa%p2+x+i+fa%p1*(y+j)
288: ENDIF
289: tonglobal=tonglobal+1
290: ENDDO
291: ENDDO
292: ENDIF
293: IF (fa%comm(0) /= 0) THEN
294: CALL DMDAGetCorners(da1,x,y,0,nx,ny,0,ierr)
295: DO j=0,ny-1
296: DO i=0,nx-1
297: tonatural(tonglobal)=fa%p1*fa%r2g+x+i+fa%p1*(y+j)
298: tonglobal=tonglobal+1
299: ENDDO
300: ENDDO
301: ENDIF
302: !-----------------------------------------------------------------------
303: ! diagnose tonatural.
304: !-----------------------------------------------------------------------
305: IF (diagnose)THEN
306: WRITE(*,'(a,i3,a)')"tonglobal = ",tonglobal,", tonatural = "
307: WRITE(*,'(2i4)')(i,tonatural(i),i=0,tonglobal-1)
308: ENDIF
309: !-----------------------------------------------------------------------
310: ! create application ordering and deallocate tonatural.
311: !-----------------------------------------------------------------------
312: CALL AOCreateBasic(PETSC_COMM_WORLD,tonglobal,tonatural, &
313: & PETSC_NULL_INTEGER,toao,ierr)
314: DEALLOCATE(tonatural)
315: !-----------------------------------------------------------------------
316: ! count the number of unknowns owned on each processor and determine
317: ! the starting point of each processors ownership for global vector
318: ! without redundancy.
319: !-----------------------------------------------------------------------
320: fromnglobal=0
321: fa%offg(1)=0
322: offt(1)=0
323: IF (fa%comm(1) /= 0)THEN
324: CALL DMDAGetCorners(da2,x,y,0,nx,ny,0,ierr)
325: offt(2)=nx*ny
326: IF (y+ny == fa%r2g)ny=ny-1
327: fromnglobal=fromnglobal+nx*ny
328: fa%offg(2)=fromnglobal
329: ELSE
330: offt(2)=0
331: fa%offg(2)=0
332: ENDIF
333: IF (fa%comm(2) /= 0)THEN
334: CALL DMDAGetCorners(da3,x,y,0,nx,ny,0,ierr)
335: offt(0)=offt(2)+nx*ny
336: IF (y+ny == fa%r2g)ny=ny-1
337: fromnglobal=fromnglobal+nx*ny
338: fa%offg(0)=fromnglobal
339: ELSE
340: offt(0)=offt(2)
341: fa%offg(0)=fromnglobal
342: ENDIF
343: IF (fa%comm(0) /= 0)THEN
344: CALL DMDAGetCorners(da1,x,y,0,nx,ny,0,ierr)
345: IF (y == 0)ny=ny-1
346: fromnglobal=fromnglobal+nx*ny
347: ENDIF
348: CALL MPI_Scan(fromnglobal,globalrstart,1,MPI_INTEGER, &
349: & MPI_SUM,PETSC_COMM_WORLD,ierr)
350: globalrstart=globalrstart-fromnglobal
351: WRITE(*,'("[",i1,"]",a,i3)') &
352: & rank," Number of unknowns owned ",fromnglobal
353: !-----------------------------------------------------------------------
354: ! get fromnatural number for each node.
355: !-----------------------------------------------------------------------
356: ALLOCATE(fromnatural(0:fromnglobal))
357: fromnglobal=0
358: IF (fa%comm(1) /= 0)THEN
359: CALL DMDAGetCorners(da2,x,y,0,nx,ny,0,ierr)
360: IF (y+ny == fa%r2g)ny=ny-1
361: fa%xg(1)=x
362: fa%yg(1)=y
363: fa%mg(1)=nx
364: fa%ng(1)=ny
365: CALL DMDAGetGhostCorners(da2,fa%xl(1),fa%yl(1),0,fa%ml(1), &
366: & fa%nl(1),0,ierr)
367: DO j=0,ny-1
368: DO i=0,nx-1
369: fromnatural(fromnglobal) &
370: & =(fa%p1-fa%p2)/2+x+i+fa%p1*(y+j)
371: fromnglobal=fromnglobal+1
372: ENDDO
373: ENDDO
374: ENDIF
375: IF (fa%comm(2) /= 0)THEN
376: CALL DMDAGetCorners(da3,x,y,0,nx,ny,0,ierr)
377: IF (y+ny == fa%r2g)ny=ny-1
378: fa%xg(2)=x
379: fa%yg(2)=y
380: fa%mg(2)=nx
381: fa%ng(2)=ny
382: CALL DMDAGetGhostCorners(da3,fa%xl(2),fa%yl(2),0,fa%ml(2), &
383: & fa%nl(2),0,ierr)
384: DO j=0,ny-1
385: DO i=0,nx-1
386: IF (x+i < (fa%p1-fa%p2)/2)THEN
387: fromnatural(fromnglobal)=x+i+fa%p1*(y+j)
388: ELSE
389: fromnatural(fromnglobal)=fa%p2+x+i+fa%p1*(y+j)
390: ENDIF
391: fromnglobal=fromnglobal+1
392: ENDDO
393: ENDDO
394: ENDIF
395: IF (fa%comm(0) /= 0)THEN
396: CALL DMDAGetCorners(da1,x,y,0,nx,ny,0,ierr)
397: IF (y == 0)THEN
398: ny=ny-1
399: ELSE
400: y=y-1
401: ENDIF
402: fa%xg(0)=x
403: fa%yg(0)=y
404: fa%mg(0)=nx
405: fa%ng(0)=ny
406: CALL DMDAGetGhostCorners(da1,fa%xl(0),fa%yl(0),0,fa%ml(0), &
407: & fa%nl(0),0,ierr)
408: DO j=0,ny-1
409: DO i=0,nx-1
410: fromnatural(fromnglobal)=fa%p1*fa%r2+x+i+fa%p1*(y+j)
411: fromnglobal=fromnglobal+1
412: ENDDO
413: ENDDO
414: ENDIF
415: !-----------------------------------------------------------------------
416: ! diagnose fromnatural.
417: !-----------------------------------------------------------------------
418: IF (diagnose)THEN
419: WRITE(*,'(a,i3,a)')"fromnglobal = ",fromnglobal &
420: & ,", fromnatural = "
421: WRITE(*,'(2i4)')(i,fromnatural(i),i=0,fromnglobal-1)
422: ENDIF
423: !-----------------------------------------------------------------------
424: ! create application ordering and deallocate fromnatural.
425: !-----------------------------------------------------------------------
426: CALL AOCreateBasic(PETSC_COMM_WORLD,fromnglobal,fromnatural, &
427: & PETSC_NULL_INTEGER,globalao,ierr)
428: DEALLOCATE(fromnatural)
429: !-----------------------------------------------------------------------
430: ! set up scatter that updates 1 from 2 and 3 and 3 and 2 from 1
431: !-----------------------------------------------------------------------
432: ALLOCATE(to(0:tonglobal),from(0:tonglobal))
433: nscat=0
434: IF (fa%comm(1) /= 0)THEN
435: CALL DMDAGetCorners(da2,x,y,0,nx,ny,0,ierr)
436: DO j=0,ny-1
437: DO i=0,nx-1
438: to(nscat)=(fa%p1-fa%p2)/2+x+i+fa%p1*(y+j)
439: from(nscat)=to(nscat)
440: nscat=nscat+1
441: ENDDO
442: ENDDO
443: ENDIF
444: IF (fa%comm(2) /= 0)THEN
445: CALL DMDAGetCorners(da3,x,y,0,nx,ny,0,ierr)
446: DO j=0,ny-1
447: DO i=0,nx-1
448: IF (x+i < (fa%p1-fa%p2)/2)THEN
449: to(nscat)=x+i+fa%p1*(y+j)
450: ELSE
451: to(nscat)=fa%p2+x+i+fa%p1*(y+j)
452: ENDIF
453: from(nscat)=to(nscat)
454: nscat=nscat+1
455: ENDDO
456: ENDDO
457: ENDIF
458: IF (fa%comm(0) /= 0)THEN
459: CALL DMDAGetCorners(da1,x,y,0,nx,ny,0,ierr)
460: DO j=0,ny-1
461: DO i=0,nx-1
462: to(nscat)=fa%p1*fa%r2g+x+i+fa%p1*(y+j)
463: from(nscat)=fa%p1*(fa%r2-1)+x+i+fa%p1*(y+j)
464: nscat=nscat+1
465: ENDDO
466: ENDDO
467: ENDIF
468: !-----------------------------------------------------------------------
469: ! diagnose to and from.
470: !-----------------------------------------------------------------------
471: IF (diagnose)THEN
472: WRITE(*,'(a,i3,a)')"nscat = ",nscat,", to, from = "
473: WRITE(*,'(3i4)')(i,to(i),from(i),i=0,nscat-1)
474: ENDIF
475: !-----------------------------------------------------------------------
476: ! create vecscatter.
477: !-----------------------------------------------------------------------
478: CALL AOApplicationToPetsc(toao,nscat,to,ierr)
479: CALL AOApplicationToPetsc(globalao,nscat,from,ierr)
480: CALL ISCreateGeneral(PETSC_COMM_WORLD,nscat,to,PETSC_COPY_VALUES, &
481: & tois,ierr)
482: CALL ISCreateGeneral(PETSC_COMM_WORLD,nscat,from,PETSC_COPY_VALUES &
483: & ,globalis,ierr)
484: CALL VecCreateMPI(PETSC_COMM_WORLD,tonglobal,PETSC_DETERMINE, &
485: & tovec,ierr)
486: CALL VecCreateMPI(PETSC_COMM_WORLD,fromnglobal,PETSC_DETERMINE, &
487: & globalvec,ierr)
488: CALL VecScatterCreate(globalvec,globalis,tovec,tois,vscat,ierr)
489: !-----------------------------------------------------------------------
490: ! clean up.
491: !-----------------------------------------------------------------------
492: CALL ISDestroy(tois,ierr)
493: CALL ISDestroy(globalis,ierr)
494: CALL AODestroy(globalao,ierr)
495: CALL AODestroy(toao,ierr)
496: DEALLOCATE(to,from)
497: !-----------------------------------------------------------------------
498: ! fill up global vector without redundant values with PETSc global
499: ! numbering
500: !-----------------------------------------------------------------------
501: CALL VecGetArray(globalvec,globalarray,ig,ierr)
502: k=ig
503: IF (diagnose)WRITE(*,'(a,i3,a)') &
504: & "fromnglobal = ",fromnglobal,", globalarray = "
505: DO i=0,fromnglobal-1
506: k=k+1
507: globalarray(k)=globalrstart+i
508: IF (diagnose)WRITE(*,'(i4,f11.3)')i,globalarray(k)
509: ENDDO
510: CALL VecRestoreArray(globalvec,globalarray,ig,ierr)
511: !-----------------------------------------------------------------------
512: ! scatter PETSc global indices to redundant valued array.
513: !-----------------------------------------------------------------------
514: CALL VecScatterBegin(vscat,globalvec,tovec,INSERT_VALUES, &
515: & SCATTER_FORWARD,ierr)
516: CALL VecScatterEnd(vscat,globalvec,tovec,INSERT_VALUES, &
517: & SCATTER_FORWARD,ierr)
518: !-----------------------------------------------------------------------
519: ! create local vector as concatenation of local vectors
520: !-----------------------------------------------------------------------
521: nlocal=0
522: cntl1=0
523: cntl2=0
524: cntl3=0
525: IF (fa%comm(1) /= 0)THEN
526: CALL VecGetSize(vl2,cntl2,ierr)
527: nlocal=nlocal+cntl2
528: ENDIF
529: IF (fa%comm(2) /= 0)THEN
530: CALL VecGetSize(vl3,cntl3,ierr)
531: nlocal=nlocal+cntl3
532: ENDIF
533: IF (fa%comm(0) /= 0)THEN
534: CALL VecGetSize(vl1,cntl1,ierr)
535: nlocal=nlocal+cntl1
536: ENDIF
537: fa%offl(0)=cntl2+cntl3
538: fa%offl(1)=0
539: fa%offl(2)=cntl2
540: CALL VecCreateSeq(PETSC_COMM_SELF,nlocal,localvec,ierr)
541: !-----------------------------------------------------------------------
542: ! cheat so that vl1, vl2, vl3 shared array memory with localvec.
543: !-----------------------------------------------------------------------
544: CALL VecGetArray(localvec,localarray,il,ierr)
545: CALL VecGetArray(tovec,toarray,it,ierr)
546: IF (fa%comm(1) /= 0)THEN
547: CALL VecPlaceArray(vl2,localarray(il+1+fa%offl(1)),ierr)
548: CALL VecPlaceArray(vg2,toarray(it+1+offt(1)),ierr)
549: CALL DMGlobalToLocalBegin(da2,vg2,INSERT_VALUES,vl2,ierr)
550: CALL DMGlobalToLocalEnd(da2,vg2,INSERT_VALUES,vl2,ierr)
551: CALL DMRestoreGlobalVector(da2,vg2,ierr)
552: ENDIF
553: IF (fa%comm(2) /= 0)THEN
554: CALL VecPlaceArray(vl3,localarray(il+1+fa%offl(2)),ierr)
555: CALL VecPlaceArray(vg3,toarray(it+1+offt(2)),ierr)
556: CALL DMGlobalToLocalBegin(da3,vg3,INSERT_VALUES,vl3,ierr)
557: CALL DMGlobalToLocalEnd(da3,vg3,INSERT_VALUES,vl3,ierr)
558: CALL DMRestoreGlobalVector(da3,vg3,ierr)
559: ENDIF
560: IF (fa%comm(0) /= 0)THEN
561: CALL VecPlaceArray(vl1,localarray(il+1+fa%offl(0)),ierr)
562: CALL VecPlaceArray(vg1,toarray(it+1+offt(0)),ierr)
563: CALL DMGlobalToLocalBegin(da1,vg1,INSERT_VALUES,vl1,ierr)
564: CALL DMGlobalToLocalEnd(da1,vg1,INSERT_VALUES,vl1,ierr)
565: CALL DMRestoreGlobalVector(da1,vg1,ierr)
566: ENDIF
567: CALL VecRestoreArray(localvec,localarray,il,ierr)
568: CALL VecRestoreArray(tovec,toarray,it,ierr)
569: !-----------------------------------------------------------------------
570: ! clean up.
571: !-----------------------------------------------------------------------
572: CALL VecScatterDestroy(vscat,ierr)
573: CALL VecDestroy(tovec,ierr)
574: !-----------------------------------------------------------------------
575: ! compute index set.
576: !-----------------------------------------------------------------------
577: ALLOCATE(indices(0:nlocal-1))
578: CALL VecGetArray(localvec,localarray,il,ierr)
579: k=il
580: IF (diagnose)WRITE(*,'(a,i3,a3,i4,a)') &
581: & "nlocal = ", nlocal,", offl = ",fa%offl,", indices = "
582: DO i=0,nlocal-1
583: k=k+1
584: indices(i)= int(localarray(k))
585: IF (diagnose)WRITE(*,'(2i4)')i,indices(i)
586: ENDDO
587: CALL VecRestoreArray(localvec,localarray,il,ierr)
588: CALL ISCreateBlock(PETSC_COMM_WORLD,2,nlocal,indices, &
589: & PETSC_COPY_VALUES,is,ierr)
590: DEALLOCATE(indices)
591: !-----------------------------------------------------------------------
592: ! create local and global vectors.
593: !-----------------------------------------------------------------------
594: CALL VecCreateSeq(PETSC_COMM_SELF,2*nlocal,fa%l,ierr)
595: CALL VecCreateMPI(PETSC_COMM_WORLD,2*fromnglobal,PETSC_DETERMINE, &
596: & fa%g,ierr)
597: !-----------------------------------------------------------------------
598: ! create final scatter that goes directly from globalvec to localvec
599: ! this is the one to be used in the application code.
600: !-----------------------------------------------------------------------
601: CALL VecScatterCreate(fa%g,is,fa%l,PETSC_NULL_OBJECT, &
602: & fa%vscat,ierr)
603: !-----------------------------------------------------------------------
604: ! clean up.
605: !-----------------------------------------------------------------------
606: CALL ISDestroy(is,ierr)
607: CALL VecDestroy(globalvec,ierr)
608: CALL VecDestroy(localvec,ierr)
609: IF (fa%comm(0) /= 0)THEN
610: CALL DMRestoreLocalVector(da1,vl1,ierr)
611: CALL DMDestroy(da1,ierr)
612: ENDIF
613: IF (fa%comm(1) /= 0)THEN
614: CALL DMRestoreLocalVector(da2,vl2,ierr)
615: CALL DMDestroy(da2,ierr)
616: ENDIF
617: IF (fa%comm(2) /= 0)THEN
618: CALL DMRestoreLocalVector(da3,vl3,ierr)
619: CALL DMDestroy(da3,ierr)
620: ENDIF
621: !-----------------------------------------------------------------------
622: ! terminate.
623: !-----------------------------------------------------------------------
624: RETURN
625: END SUBROUTINE barry_create_fa
626: !-----------------------------------------------------------------------
627: ! subprogram 7. barry_draw_patch.
628: ! crude graphics to test that the ghost points are properly updated.
629: !-----------------------------------------------------------------------
630: !-----------------------------------------------------------------------
631: ! declarations.
632: !-----------------------------------------------------------------------
633: SUBROUTINE barry_draw_patch(draw,patch,ierr)
635: PetscDraw, INTENT(IN) :: draw
636: TYPE(patch_type), DIMENSION(0:2), INTENT(IN) :: patch
637: INTEGER, INTENT(OUT) :: ierr
639: INTEGER :: ix,iy,j
640: PetscReal, DIMENSION(4) :: x,y
641: !-----------------------------------------------------------------------
642: ! draw it.
643: !-----------------------------------------------------------------------
644: DO j=0,2
645: DO iy=1,patch(j)%my
646: DO ix=1,patch(j)%mx
647: x(1)=patch(j)%xy(1,ix-1,iy-1)
648: y(1)=patch(j)%xy(2,ix-1,iy-1)
649: x(2)=patch(j)%xy(1,ix,iy-1)
650: y(2)=patch(j)%xy(2,ix,iy-1)
651: x(3)=patch(j)%xy(1,ix,iy)
652: y(3)=patch(j)%xy(2,ix,iy)
653: x(4)=patch(j)%xy(1,ix-1,iy)
654: y(4)=patch(j)%xy(2,ix-1,iy)
655: CALL PetscDrawLine(draw,x(1),y(1),x(2),y(2), &
656: & PETSC_DRAW_BLACK,ierr)
657: CALL PetscDrawLine(draw,x(2),y(2),x(3),y(3), &
658: & PETSC_DRAW_BLACK,ierr)
659: CALL PetscDrawLine(draw,x(3),y(3),x(4),y(4), &
660: & PETSC_DRAW_BLACK,ierr)
661: CALL PetscDrawLine(draw,x(4),y(4),x(1),y(1), &
662: & PETSC_DRAW_BLACK,ierr)
663: ENDDO
664: ENDDO
665: ENDDO
666: !-----------------------------------------------------------------------
667: ! terminate.
668: !-----------------------------------------------------------------------
669: ierr=0
670: RETURN
671: END SUBROUTINE barry_draw_patch
672: !-----------------------------------------------------------------------
673: ! subprogram 8. barry_draw_fa.
674: ! deallocates local array.
675: !-----------------------------------------------------------------------
676: !-----------------------------------------------------------------------
677: ! declarations.
678: !-----------------------------------------------------------------------
679: SUBROUTINE barry_draw_fa(fa,v)
681: TYPE(fa_type), INTENT(IN) :: fa
682: Vec, INTENT(IN) :: v
684: INTEGER :: iv,vn,ln,j,offset
685: PetscReal, DIMENSION(1) :: va
686: TYPE(patch_type), DIMENSION(0:2) :: patch
687: PetscDraw :: draw
688: PetscReal :: xmin,xmax,ymin,ymax
689: PetscReal :: ymint,xmint,ymaxt,xmaxt
690: !-----------------------------------------------------------------------
691: ! set extrema.
692: !-----------------------------------------------------------------------
693: xmin=HUGE(xmin)
694: xmax=-HUGE(xmax)
695: ymin=HUGE(ymin)
696: ymax=-HUGE(ymax)
697: xmint=HUGE(xmint)
698: xmaxt=-HUGE(xmaxt)
699: ymint=HUGE(ymint)
700: ymaxt=-HUGE(ymaxt)
701: !-----------------------------------------------------------------------
702: ! get arrays and sizes.
703: !-----------------------------------------------------------------------
704: CALL VecGetArray(v,va,iv,ierr)
705: CALL VecGetSize(v,vn,ierr)
706: CALL VecGetSize(fa%l,ln,ierr)
707: !-----------------------------------------------------------------------
708: ! fill arrays.
709: !-----------------------------------------------------------------------
710: DO j=0,2
711: IF (vn == ln)THEN
712: offset=iv+2*fa%offl(j)
713: patch(j)%mx=fa%ml(j)-1
714: patch(j)%my=fa%nl(j)-1
715: ELSE
716: offset=iv+2*fa%offg(j)
717: patch(j)%mx=fa%mg(j)-1
718: patch(j)%my=fa%ng(j)-1
719: ENDIF
720: ALLOCATE(patch(j)%xy(2,0:patch(j)%mx,0:patch(j)%my))
721: patch(j)%xy=RESHAPE(va(offset+1:offset+SIZE(patch(j)%xy)), &
722: & (/2,patch(j)%mx+1,patch(j)%my+1/))
723: ENDDO
724: !-----------------------------------------------------------------------
725: ! compute extrema over processor.
726: !-----------------------------------------------------------------------
727: DO j=0,2
728: xmin=MIN(xmin,MINVAL(patch(j)%xy(1,:,:)))
729: xmax=MAX(xmax,MAXVAL(patch(j)%xy(1,:,:)))
730: ymin=MIN(ymin,MINVAL(patch(j)%xy(2,:,:)))
731: ymax=MAX(ymax,MAXVAL(patch(j)%xy(2,:,:)))
732: ENDDO
733: !-----------------------------------------------------------------------
734: ! compute global extrema.
735: !-----------------------------------------------------------------------
736: CALL MPI_Allreduce(xmin,xmint,1,MPI_DOUBLE_PRECISION,MPI_MIN, &
737: & PETSC_COMM_WORLD,ierr)
738: CALL MPI_Allreduce(xmax,xmaxt,1,MPI_DOUBLE_PRECISION,MPI_MAX, &
739: & PETSC_COMM_WORLD,ierr)
740: CALL MPI_Allreduce(ymin,ymint,1,MPI_DOUBLE_PRECISION,MPI_MIN, &
741: & PETSC_COMM_WORLD,ierr)
742: CALL MPI_Allreduce(ymax,ymaxt,1,MPI_DOUBLE_PRECISION,MPI_MAX, &
743: & PETSC_COMM_WORLD,ierr)
744: !-----------------------------------------------------------------------
745: ! set margins.
746: !-----------------------------------------------------------------------
747: xmin=xmint-.2*(xmaxt-xmint)
748: xmax=xmaxt+.2*(xmaxt-xmint)
749: ymin=ymint-.2*(ymaxt-ymint)
750: ymax=ymaxt+.2*(ymaxt-ymint)
751: !-----------------------------------------------------------------------
752: ! draw it.
753: !-----------------------------------------------------------------------
754: #if defined(PETSC_HAVE_X) || defined(PETSC_HAVE_OPENGL)
755: CALL PetscDrawCreate(PETSC_COMM_WORLD,0,"meshes",PETSC_DECIDE, &
756: & PETSC_DECIDE,700,700,draw,ierr)
757: CALL PetscDrawSetFromOptions(draw,ierr)
758: CALL PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax,ierr)
759: CALL PetscDrawZoom(draw,barry_draw_patch,patch,ierr)
760: #endif
761: !-----------------------------------------------------------------------
762: ! clean up.
763: !-----------------------------------------------------------------------
764: CALL VecRestoreArray(v,va,iv,ierr)
765: #if defined(PETSC_HAVE_X)
766: CALL PetscDrawDestroy(draw,ierr)
767: #endif
768: !-----------------------------------------------------------------------
769: ! terminate.
770: !-----------------------------------------------------------------------
771: RETURN
772: END SUBROUTINE barry_draw_fa
773: !-----------------------------------------------------------------------
774: ! subprogram 9. barry_map_region3.
775: ! fills region 3 coordinates.
776: !-----------------------------------------------------------------------
777: !-----------------------------------------------------------------------
778: ! declarations.
779: !-----------------------------------------------------------------------
780: SUBROUTINE barry_map_region3(fa,g)
782: TYPE(fa_type), INTENT(INOUT) :: fa
783: Vec, INTENT(INOUT) :: g
785: INTEGER :: i,j,k,x,y,m,n,ig
786: PetscReal :: r0=1,theta0=pi/6,r,theta
787: PetscReal :: dr,dt
788: PetscReal, DIMENSION(1) :: ga
789: !-----------------------------------------------------------------------
790: ! format statements.
791: !-----------------------------------------------------------------------
792: 10 FORMAT(/5x,"j",5x,"i",6x,"r",8x,"theta",8x,"x",10x,"y"/)
793: 20 FORMAT(2i6,4f11.3)
794: !-----------------------------------------------------------------------
795: ! preliminary computations.
796: !-----------------------------------------------------------------------
797: dr=r0/(fa%r2-1)
798: dt=twopi/(3*(fa%p1-fa%p2-1))
799: CALL barry_get_global_corners(fa,2,x,y,m,n)
800: !-----------------------------------------------------------------------
801: ! fill array.
802: !-----------------------------------------------------------------------
803: CALL VecGetArray(g,ga,ig,ierr)
804: k=ig+2*fa%offg(2)
805: IF (diagnose)THEN
806: WRITE(*,'(a)')"region 3"
807: WRITE(*,10)
808: ENDIF
809: DO j=y,y+n-1
810: r=r0+j*dr
811: DO i=x,x+m-1
812: theta=theta0+i*dt
813: ga(k+1)=r*COS(theta)
814: ga(k+2)=r*SIN(theta)-4*r0
815: IF (diagnose)WRITE(*,20)j,i,r,theta,ga(k+1),ga(k+2)
816: k=k+2
817: ENDDO
818: IF (diagnose)WRITE(*,10)
819: ENDDO
820: CALL VecRestoreArray(g,ga,ig,ierr)
821: !-----------------------------------------------------------------------
822: ! terminate.
823: !-----------------------------------------------------------------------
824: RETURN
825: END SUBROUTINE barry_map_region3
826: !-----------------------------------------------------------------------
827: ! subprogram 10. barry_map_region2.
828: ! fills region 2 coordinates.
829: !-----------------------------------------------------------------------
830: !-----------------------------------------------------------------------
831: ! declarations.
832: !-----------------------------------------------------------------------
833: SUBROUTINE barry_map_region2(fa,g)
835: TYPE(fa_type), INTENT(INOUT) :: fa
836: Vec, INTENT(INOUT) :: g
838: INTEGER :: i,j,k,x,y,m,n,ig
839: PetscReal :: r0=1,theta0=-pi/2,r,theta
840: PetscReal :: dr,dt
841: PetscReal, DIMENSION(1) :: ga
842: !-----------------------------------------------------------------------
843: ! format statements.
844: !-----------------------------------------------------------------------
845: 10 FORMAT(/5x,"j",5x,"i",6x,"r",8x,"theta",8x,"x",10x,"y"/)
846: 20 FORMAT(2i6,4f11.3)
847: !-----------------------------------------------------------------------
848: ! preliminary computations.
849: !-----------------------------------------------------------------------
850: dr=r0/(fa%r2-1)
851: dt=twopi/fa%p2
852: CALL barry_get_global_corners(fa,1,x,y,m,n)
853: !-----------------------------------------------------------------------
854: ! fill array.
855: !-----------------------------------------------------------------------
856: CALL VecGetArray(g,ga,ig,ierr)
857: k=ig+2*fa%offg(1)
858: IF (diagnose)THEN
859: WRITE(*,'(a)')"region 2"
860: WRITE(*,10)
861: ENDIF
862: DO j=y,y+n-1
863: r=r0+j*dr
864: DO i=x,x+m-1
865: theta=theta0+i*dt
866: ga(k+1)=r*COS(theta)
867: ga(k+2)=r*SIN(theta)
868: IF (diagnose)WRITE(*,20)j,i,r,theta,ga(k+1),ga(k+2)
869: k=k+2
870: ENDDO
871: IF (diagnose)WRITE(*,10)
872: ENDDO
873: CALL VecRestoreArray(g,ga,ig,ierr)
874: !-----------------------------------------------------------------------
875: ! terminate.
876: !-----------------------------------------------------------------------
877: RETURN
878: END SUBROUTINE barry_map_region2
879: !-----------------------------------------------------------------------
880: ! subprogram 11. barry_map_region1.
881: ! fills region 1 coordinates.
882: !-----------------------------------------------------------------------
883: !-----------------------------------------------------------------------
884: ! declarations.
885: !-----------------------------------------------------------------------
886: SUBROUTINE barry_map_region1(fa,g)
888: TYPE(fa_type), INTENT(INOUT) :: fa
889: Vec, INTENT(INOUT) :: g
891: INTEGER :: i,j,k,x,y,m,n,ig
892: PetscReal :: r0=1,r,theta,dr,dt1,dt3
893: PetscReal, DIMENSION(1) :: ga
894: !-----------------------------------------------------------------------
895: ! format statements.
896: !-----------------------------------------------------------------------
897: 10 FORMAT(/5x,"j",5x,"i",6x,"r",8x,"theta",8x,"x",10x,"y"/)
898: 20 FORMAT(2i6,4f11.3)
899: !-----------------------------------------------------------------------
900: ! preliminary computations.
901: !-----------------------------------------------------------------------
902: dr=r0/(fa%r1-1)
903: dt1=twopi/fa%p2
904: dt3=twopi/(3*(fa%p1-fa%p2-1))
905: CALL barry_get_global_corners(fa,0,x,y,m,n)
906: !-----------------------------------------------------------------------
907: ! fill array.
908: !-----------------------------------------------------------------------
909: CALL VecGetArray(g,ga,ig,ierr)
910: k=ig+2*fa%offg(0)
911: IF (diagnose)THEN
912: WRITE(*,'(a)')"region 1"
913: WRITE(*,10)
914: ENDIF
915: DO j=y,y+n-1
916: r=2*r0+j*dr
917: DO i=x,x+m-1
918: IF (i < (fa%p1-fa%p2)/2)THEN
919: theta=i*dt3
920: ga(k+1)=r*COS(theta)
921: ga(k+2)=r*SIN(theta)-4*r0
922: ELSEIF (i > fa%p2 + (fa%p1-fa%p2)/2)then
923: theta=pi+i*dt3
924: ga(k+1)=r*COS(PETSC_PI+i*dt3)
925: ga(k+2)=r*SIN(PETSC_PI+i*dt3)-4*r0
926: ELSE
927: theta=(i-(fa%p1-fa%p2)/2)*dt1-pi/2
928: ga(k+1)=r*COS(theta)
929: ga(k+2)=r*SIN(theta)
930: ENDIF
931: IF (diagnose)WRITE(*,20)j,i,r,theta,ga(k+1),ga(k+2)
932: k=k+2
933: ENDDO
934: IF (diagnose)WRITE(*,10)
935: ENDDO
936: CALL VecRestoreArray(g,ga,ig,ierr)
937: !-----------------------------------------------------------------------
938: ! terminate.
939: !-----------------------------------------------------------------------
940: RETURN
941: END SUBROUTINE barry_map_region1
942: END MODULE barry_mod
943: !-----------------------------------------------------------------------
944: ! subprogram 12. barry_main.
945: ! main program.
946: !-----------------------------------------------------------------------
947: !-----------------------------------------------------------------------
948: ! declarations.
949: !-----------------------------------------------------------------------
950: PROGRAM barry_main
951: USE barry_mod
952: IMPLICIT NONE
954: TYPE(fa_type) :: fa
955: Vec :: g,l
956: !-----------------------------------------------------------------------
957: ! initialize and create structure, and get vectors.
958: !-----------------------------------------------------------------------
959: CALL PetscInitialize(PETSC_NULL_CHARACTER,ierr)
960: CALL barry_create_fa(fa)
961: CALL barry_get_global_vector(fa,g)
962: CALL barry_get_local_vector(fa,l)
963: !-----------------------------------------------------------------------
964: ! map regions.
965: !-----------------------------------------------------------------------
966: CALL barry_map_region1(fa,g)
967: CALL barry_map_region2(fa,g)
968: CALL barry_map_region3(fa,g)
969: !-----------------------------------------------------------------------
970: ! graphic output.
971: !-----------------------------------------------------------------------
972: CALL barry_global_to_local(fa,g,l)
973: CALL barry_draw_fa(fa,g)
974: CALL barry_draw_fa(fa,l)
975: !-----------------------------------------------------------------------
976: ! clean up and finalize.
977: !-----------------------------------------------------------------------
978: CALL VecDestroy(g,ierr)
979: CALL VecDestroy(l,ierr)
980: CALL barry_destroy_fa(fa)
981: CALL PetscFinalize(ierr)
982: !-----------------------------------------------------------------------
983: ! terminate.
984: !-----------------------------------------------------------------------
985: END PROGRAM barry_main