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