GRIDMAN
grid managment library
carre.f
Go to the documentation of this file.
1 C> @file convert/carre.f
2 C> Converter for B2 (CARRE, SONNET) grid into GRIDMAN_GRID object
3 C GRIDMAN, grid managment library. Author: Vladislav Kotov, v.kotov@fz-juelich.de
4 
5 ! Copyright (c) 2017 Forschungszentrum Juelich GmbH
6 ! Vladislav Kotov
7 !
8 ! This file is part of GRIDMAN.
9 !
10 ! GRIDMAN is free software: you can redistribute it and/or modify
11 ! it under the terms of the GNU General Public License as published by
12 ! the Free Software Foundation, either version 3 of the License, or
13 ! (at your option) any later version.
14 !
15 ! GRIDMAN is distributed in the hope that it will be useful,
16 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
17 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 ! GNU General Public License for more details.
19 !
20 ! You should have received a copy of the GNU General Public License
21 ! along with GRIDMAN. If not, see <http://www.gnu.org/licenses/>.
22 
23 C> Convert B2 (CARRE, SONNET) grid into GRIDMAN_GRID grid object
24 C>
25 C> WARNING: the object GRID will be overwritten if already exists!
26 C>
27 C> 3rd index in BR, BZ is the index of corner
36 C> This numbering of corners is used BR, BZ arrays produced by
37 C> \ref GRIDMAN_CARRE_READ30, \ref GRIDMAN_CARRE_READSONNET
38 C>
39 C> Arrays NXCUT,NYCUT,NXISO,NYISO describe cuts an isolated cells as
40 C> returned by subroutine GRIDMAN_CARRE_READSONNET.
41 C> In order to use the output of GRIDMAN_CARRE_READ30
42 C> this subroutine has one more entry:
43 C>
44 C> GRIDMAN_CARRE2GRID_FORT30(GRID,LX,LY,NXL,NYL,NNCUT,
45 C> NXCUT30,NYCUT30,
46 C> NNISO,NXISO30,NYISO30,BR,BZ,IERR)
47 C>
48 C>
49 C> In the resulting grid one cell index, and 3 edge indexes.
50 C> Cell index show IX,IY indexes of the correcponding cell on the plasma (B2) grid.
51 C> 1st edge index contains markers of the grid boundaries:
52 C> -1 stands for South, -4 for North, -2 for West, -3 for East
53 C> 2nd edge index is the indices of the poloidal edges:
54 C>
55 C> _________________
56 C> | |
57 C> | o IX,IY iedge -> IX,IY
58 C> |_________________|
59 C>
60 C> 3rd edge index is the indices of the radial edges:
61 C>
62 C>
63 C> __iedge -> IX,IY _
64 C> | |
65 C> | o IX,IY |
66 C> |__________________|
67 C>
68  SUBROUTINE gridman_carre2grid(GRID,
69  s lxb,lyb,nxb,nyb,nncut,nxcut,nycut,
70  s nniso,nxiso,nyiso,br,bz,
71  s lx,ly,nxl,nyl,ierr,leir)
73  u gridman_unit,gridman_dbg,gridman_check
76  IMPLICIT NONE
77  INTRINSIC floor,mod,abs
78 
79 C> Resulting 2D grid object
80  TYPE(gridman_grid) :: GRID
81 C> First index in arrays BR, BZ for X-direction
82  INTEGER(GRIDMAN_SP),INTENT(IN) :: LXB
83 C> First index in arrays BR, BZ for Y-direction
84  INTEGER(GRIDMAN_SP),INTENT(IN) :: LYB
85 C> Last index in arrays BR, BZ for X-direction
86  INTEGER(GRIDMAN_SP),INTENT(IN) :: NXB
87 C> Last index in arrays BR, BZ for Y-direction
88  INTEGER(GRIDMAN_SP),INTENT(IN) :: NYB
89 C> Number of cuts
90  INTEGER(GRIDMAN_SP),INTENT(IN) :: NNCUT
91 C> Number of isolated surfaces
92  INTEGER(GRIDMAN_SP),INTENT(IN) :: NNISO
93  INTEGER(GRIDMAN_SP),INTENT(IN) :: NXCUT30(2,nncut) !COORDINATES OF THE CUTS, FORT30 FORMAT
94  INTEGER(GRIDMAN_SP),INTENT(IN) :: NYCUT30(2,nncut) !
95  INTEGER(GRIDMAN_SP),INTENT(IN) :: NXISO30(2,nniso) !COORDINATES OF THE ISOLATED SURFACES, ...
96  INTEGER(GRIDMAN_SP),INTENT(IN) :: NYISO30(2,nniso) !
97 C> X-indices of cuts (in GRIDMAN_CARRE_READSONNET format)
98  INTEGER(GRIDMAN_SP),INTENT(IN) :: NXCUT(nncut)
99 C> Y-indices of cuts (in GRIDMAN_CARRE_READSONNET format)
100  INTEGER(GRIDMAN_SP),INTENT(IN) :: NYCUT(nncut)
101 C> X-indices of isolated cells (in GRIDMAN_CARRE_READSONNET format)
102  INTEGER(GRIDMAN_SP),INTENT(IN) :: NXISO(nniso)
103 C> Y-indices of isolated cells (in GRIDMAN_CARRE_READSONNET format)
104  INTEGER(GRIDMAN_SP),INTENT(IN) :: NYISO(nniso)
105 C> Major radius (R) of four corners of each cell, BR(LXB:NXB,LYB:NYB,4)
106  REAL(GRIDMAN_DP) :: BR(lxb:nxb,lyb:nyb,4)
107 C> Vertical coordinate (Z) of four corners of each cell, BZ(LXB:NXB,LYB:NYB,4)
108  REAL(GRIDMAN_DP) :: BZ(lxb:nxb,lyb:nyb,4)
109 C> 1st index to be taken in arrays BR, BZ for X-direction
110  INTEGER(GRIDMAN_SP),INTENT(IN) :: LX
111 C> 1st index to be taken in arrays BR, BZ for Y-direction
112  INTEGER(GRIDMAN_SP),INTENT(IN) :: LY
113 C> Last index to be taken in arrays BR, BZ for X-direction
114  INTEGER(GRIDMAN_SP),INTENT(IN) :: NXL
115 C> Last index to be taken in arrays BR, BZ for Y-direction
116  INTEGER(GRIDMAN_SP),INTENT(IN) :: NYL
117 C> Error code
118  INTEGER,INTENT(OUT) :: IERR
119 C> Switch to produce cell index mapping in the sense of EIRENE
120 C>
121 C> If .TRUE. then inices in the index mapping of cells are shifted
122 C> in the presence if cuts in the way it is done in EIRENE
123  LOGICAL,OPTIONAL :: LEIR
124 
125  INTEGER(GRIDMAN_SP) :: NX,NY,NCELLS,NPOL,NRAD,NEDGES,NPOINTS,NBND,
126  i icell,ix,iy,ic,ip,ie,m,iedge,
127  i ipoint,icell1,icell2,ix1,ix2,iy1,iy2,
128  i ipoint1,ipoint2
129  INTEGER :: IERR0,RES
130  REAL(GRIDMAN_DP) :: X1,X2,Y1,Y2
131  LOGICAL :: LFORT30,LFREEEDGES
132 
133  lfort30=.false.
134 
135  GOTO 10
136 
137 C THIS ENTRY READS THE COORDINATES OF CUTS AND ISOLATED CELLS
138 C IN FORT.30 FORMAT, INSTEAD OF CARRE FORMAT
139  entry gridman_carre2grid_fort30(grid,lxb,lyb,nxb,nyb,nncut,
140  s nxcut30,nycut30,nniso,
141  s nxiso30,nyiso30,br,bz,
142  s lx,ly,nxl,nyl,ierr,leir)
143 
144  lfort30=.true.
145 
146  10 CONTINUE
147 
148  IF(gridman_dbg)
149  w WRITE(gridman_unit,*) "Starting GRIDMAN_CARRE2GRID"
150 
151  ierr=0
152 
153  IF(nxl.LT.lx+1.OR.nyl.LT.ly+1) THEN
154  ierr=100
155  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
156  w "wrong dimensions, LX, LY, NX, NY ",
157  w lx, ly, nx, ny
158  RETURN
159  END IF
160 
161 C CHECK CUTS
162  IF(lfort30) THEN
163  DO ic=1,nncut
164  ix1=nxcut30(1,ic)
165  ix2=nxcut30(2,ic)
166  IF(ix1.LT.lx.OR.ix2.LT.lx.OR.ix1.GT.nxl.OR.ix2.GT.nxl.OR.
167  f nycut30(2,ic).LT.nycut30(1,ic).OR.nycut30(2,ic).GT.nyl) THEN
168  ierr=100
169  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
170  o "incorrect parameters of vertical cuts"
171  WRITE(gridman_unit,*) " NXCUT1 ", nxcut30(1,:)
172  WRITE(gridman_unit,*) " NXCUT2 ", nxcut30(2,:)
173  WRITE(gridman_unit,*) " NYCUT1 ", nycut30(1,:)
174  WRITE(gridman_unit,*) " NYCUT2 ", nycut30(2,:)
175  WRITE(gridman_unit,*) " MIN(IX), MAX(IX), MIN(IY), MAX(IY) ",
176  w lx,nxl,ly,nyl
177  RETURN
178  END IF
179  END DO
180  ELSE
181  DO ic=1,nncut
182  ix1=nxcut(ic)
183  ix2=nxcut(ic)
184  IF(ix1.LT.lx-1.OR.ix2.LT.lx-1.OR.ix1.GT.nxl.OR.ix2.GT.nxl.OR.
185  f nycut(ic).LT.nycut(ic).OR.nycut(ic).GT.nyl) THEN
186  ierr=100
187  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
188  o "incorrect parameters of vertical cuts"
189  WRITE(gridman_unit,*) " NXCUT ", nxcut(1:nncut)
190  WRITE(gridman_unit,*) " NYCUT ", nycut(1:nncut)
191  WRITE(gridman_unit,*) " MIN(IX), MAX(IX), MIN(IY), MAX(IY) ",
192  w lx,nxl,ly,nyl
193  RETURN
194  END IF
195  END DO
196  END IF
197 
198  nx=nxl-lx+1
199  ny=nyl-ly+1
200  ncells=nx*ny !NUMBER OF CELLS
201  npol=nx*(ny+1) !NUMBER OF POLOIDAL EDGES
202  nrad=ny*(nx+1) !NUMBER OF RADIAL EDGES
203  nedges=npol+nrad !NUMBER OF EDGES
204  npoints=(nx+1)*(ny+1) !NUMBER OF POINTS
205  nbnd=2*nx+2*ny !TOTAL NUMBER OF EDGES ON OPEN BOUNDARIES
206 
207  m=npol
208 
209  CALL gridman_grid_allocate(grid,2,nedges,npoints,ncells,ierr,3,1)
210  IF(ierr.NE.0) GOTO 1000
211  CALL gridman_index_allocate(grid%EDGEINDEX(1),1,1_gridman_sp,ierr) !temporarary
212  IF(ierr.NE.0) THEN
213  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
214  w "cannot allocate temporar 1st edge index"
215  GOTO 2000
216  END IF
217  grid%EDGEINDEX(1)%INDEXES(0,1)=1
218  grid%EDGEINDEX(1)%INDEXES(1,1)=1
219  CALL gridman_index_allocate(grid%EDGEINDEX(2),2,npol,ierr)
220  IF(ierr.NE.0) THEN
221  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
222  w "cannot allocate index of poloidal edges"
223  GOTO 2000
224  END IF
225  CALL gridman_index_allocate(grid%EDGEINDEX(3),2,nrad,ierr)
226  IF(ierr.NE.0) THEN
227  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
228  w "cannot allocate index of radial edges"
229  GOTO 2000
230  END IF
231  CALL gridman_index_allocate(grid%CELLINDEX(1),2,ncells,ierr)
232  IF(ierr.NE.0) THEN
233  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
234  w "cannot allocate index of cells"
235  GOTO 2000
236  END IF
237 C
238 C 3---------2 |----2----|
239 C | | | GRIDMAN |
240 C ^ | B2 CELL | 3 1
241 C YI | | | CELL |
242 C 4---------1 |----4----|
243 C XI >
244 C
245 
246 C
247 C EXAMPLE OF A GRID
248 C
249 C NX=4, NY=3, M=NX*(NY+1)=16
250 C
251 C 16--13---17--14---18--15---19--16---20
252 C | | | | |
253 C 27 9 28 10 29 11 30 12 31
254 C | | | | |
255 C 11--9----12--10---13--11---14--12---15
256 C | | | | |
257 C 22 5 23 6 24 7 25 8 26
258 C | | | | |
259 C 6---5----7---6----8---7----8---8----10
260 C | | | | |
261 C 17 1 18 2 19 3 20 4 21
262 C | | | | |
263 C 1---1----2---2----3---3----4---4----5
264 C
265 
266 C
267 C IX,IY IN THE CYCLES ARE SHIFTED COORDINATES: WHICH STARTS FROM LX,LY
268 C THEY ARE RELATED TO THE COORDINATES WHICH START FROM 1 AS
269 C IX_LX = IX_1+LX-1 , IY_LY = IY_1+LY-1
270 C
271 
272 
273 C POLOIDAL (X) EDGES
274 C
275 C | |
276 C | IEDGE |
277 C | |
278 C IEDGE+IY-1 o------ IEDGE ----- o IEDGE+IY
279 C | |
280 C | IEDGE-NX |
281 C | |
282 C
283  iedge=0
284  DO ix=lx,nxl
285  iedge=iedge+1
286  grid%CELLS(1,iedge)=iedge
287  grid%CELLS(2,iedge)=0
288  grid%POINTS(1,iedge)=iedge
289  grid%POINTS(2,iedge)=grid%POINTS(1,iedge)+1
290  grid%EDGEINDEX(2)%INDEXES(0,iedge)=iedge
291  grid%EDGEINDEX(2)%INDEXES(1,iedge)=ix
292  grid%EDGEINDEX(2)%INDEXES(2,iedge)=ly-1
293  END DO
294  DO iy=ly+1,nyl
295  DO ix=lx,nxl
296  iedge=iedge+1
297  grid%CELLS(1,iedge)=iedge-nx
298  grid%CELLS(2,iedge)=iedge
299  grid%POINTS(1,iedge)=iedge+iy-ly
300  grid%POINTS(2,iedge)=grid%POINTS(1,iedge)+1
301  grid%EDGEINDEX(2)%INDEXES(0,iedge)=iedge
302  grid%EDGEINDEX(2)%INDEXES(1,iedge)=ix
303  grid%EDGEINDEX(2)%INDEXES(2,iedge)=iy-1
304  END DO
305  END DO
306  DO ix=lx,nxl
307  iedge=iedge+1
308  grid%CELLS(1,iedge)=iedge-nx
309  grid%CELLS(2,iedge)=0
310  grid%POINTS(1,iedge)=iedge+ny
311  grid%POINTS(2,iedge)=grid%POINTS(1,iedge)+1
312  grid%EDGEINDEX(2)%INDEXES(0,iedge)=iedge
313  grid%EDGEINDEX(2)%INDEXES(1,iedge)=ix
314  grid%EDGEINDEX(2)%INDEXES(2,iedge)=nyl
315  END DO
316  IF(iedge.NE.m) THEN
317  ierr=400
318  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
319  o "specification of CELLS is incorrect"
320  WRITE(gridman_unit,*) " IEDGE, M, NY, NX, ",
321  o iedge, m, ny, nx
322  GOTO 2000
323  END IF
324 
325 C RADIAL (Y) EDGES
326 C
327 C IEDGE-M+NX+1
328 C -----------------o------------------
329 C |
330 C IEDGE-M-IY IEDGE IEDGE-M-IY+1
331 C |
332 C -----------------o------------------
333 C IEDGE-M
334 C
335  ie=0 !IE=IEDGE-NPOL
336  DO iy=ly,nyl
337  iedge=iedge+1
338  grid%CELLS(1,iedge)=iedge-m-iy+ly
339  grid%CELLS(2,iedge)=0
340  grid%POINTS(1,iedge)=iedge-m
341  grid%POINTS(2,iedge)=grid%POINTS(1,iedge)+nx+1
342  ie=ie+1
343  grid%EDGEINDEX(3)%INDEXES(0,ie)=iedge
344  grid%EDGEINDEX(3)%INDEXES(1,ie)=lx-1
345  grid%EDGEINDEX(3)%INDEXES(2,ie)=iy
346  DO ix=lx+1,nxl
347  iedge=iedge+1
348  grid%CELLS(1,iedge)=iedge-m-(iy-ly+1)
349  grid%CELLS(2,iedge)=iedge-m-(iy-ly+1)+1
350  grid%POINTS(1,iedge)=iedge-m
351  grid%POINTS(2,iedge)=grid%POINTS(1,iedge)+nx+1
352  ie=ie+1
353  grid%EDGEINDEX(3)%INDEXES(0,ie)=iedge
354  grid%EDGEINDEX(3)%INDEXES(1,ie)=ix-1
355  grid%EDGEINDEX(3)%INDEXES(2,ie)=iy
356  END DO
357  iedge=iedge+1
358  grid%CELLS(1,iedge)=iedge-m-(iy-ly+1)
359  grid%CELLS(2,iedge)=0
360  grid%POINTS(1,iedge)=iedge-m
361  grid%POINTS(2,iedge)=grid%POINTS(1,iedge)+nx+1
362  ie=ie+1
363  grid%EDGEINDEX(3)%INDEXES(0,ie)=iedge
364  grid%EDGEINDEX(3)%INDEXES(1,ie)=nxl
365  grid%EDGEINDEX(3)%INDEXES(2,ie)=iy
366  END DO
367  IF(iedge.NE.nedges) THEN
368  ierr=400
369  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
370  o "specification of CELLS is incorrect"
371  WRITE(gridman_unit,*) " IEDGE, NEDGES, NY, NX ",
372  o iedge, nedges, ny, nx
373  GOTO 2000
374  END IF
375 C
376 C COORDINATES
377 C
378 C IY
379 C 9___10___11___12
380 C |\ |\ |\ /|
381 C 2 | 32 | 32 | 32 |
382 C | 41 | 41 | 41 |
383 C |/___|/__ |/__\|
384 C 5 6 7 8
385 C | | | |
386 C 1 | 32 | 32 | 32 |
387 C | 41 | 41 | 41 |
388 C |/___|/___|/__\|
389 C 1 2 3 4
390 C IX 1 2 3
391 C
392 C
393  ipoint=0
394  DO iy=ly,nyl
395  DO ix=lx,nxl
396  ipoint=ipoint+1
397  grid%X(1,ipoint)=br(ix,iy,4)
398  grid%X(2,ipoint)=bz(ix,iy,4)
399  END DO
400  ipoint=ipoint+1
401  grid%X(1,ipoint)=br(nxl,iy,1)
402  grid%X(2,ipoint)=bz(nxl,iy,1)
403  END DO
404  DO ix=lx,nxl
405  ipoint=ipoint+1
406  grid%X(1,ipoint)=br(ix,nyl,3)
407  grid%X(2,ipoint)=bz(ix,nyl,3)
408  END DO
409  ipoint=ipoint+1
410  grid%X(1,ipoint)=br(nxl,nyl,2)
411  grid%X(2,ipoint)=bz(nxl,nyl,2)
412  IF(ipoint.NE.npoints) THEN
413  ierr=400
414  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
415  o "specification of points is incorrect"
416  WRITE(gridman_unit,*) " IPOINT, NPOINTS, NY, NX, ",
417  o ipoint, npoints, ny, nx
418  GOTO 2000
419  END IF
420 
421 C INDEXES OF CELLS
422  DO iy=ly,nyl
423  DO ix=lx,nxl
424  icell=ix-lx+1+nx*(iy-ly)
425  grid%CELLINDEX(1)%INDEXES(0,icell)=icell
426  grid%CELLINDEX(1)%INDEXES(1,icell)=ix
427  grid%CELLINDEX(1)%INDEXES(2,icell)=iy
428  END DO
429  END DO
430 
431  m=nx*(ny+1) !INDEX SHIFT FOR VERTICAL EDGES
432 C
433 C ADD VERTICAL CUTS
434 C
435 
436 C TRANSFORMATION OF INDICES
437 C
438 C ICELL=IX+NX*(IY-1)
439 C IEDGE4=ICELL
440 C IEDGE2=ICELL+NX
441 C IEDGE1=ICELL+M+IY-LY+1
442 C IEDGE3=ICELL+M+IY-LY
443 C
444  lfreeedges=.false.
445  DO ic=1,nncut
446  IF(lfort30) THEN
447  ix1=nxcut30(1,ic)
448  ix2=nxcut30(2,ic)
449  iy1=max(nycut30(1,ic),ly)
450  iy2=nycut30(2,ic)
451  ELSE
452  ip=nncut-ic+1
453  ix1=nxcut(ic)
454  ix2=nxcut(ip)+1
455  iy1=ly
456  iy2=nycut(ic)
457  END IF
458 C SPECIAL TREATMENT OF THE CASE W/O DIVERTOR LEGS
459  IF(ix1.LT.lx) THEN
460  DO iy=iy1,iy2
461  icell2=ix2-1-lx+1+nx*(iy-ly) !(IX2-1,IY)
462  iedge=icell2+iy-ly+m+1 !SIDE 1 OF (IX2-1,IY)
463  grid%CELLS(1,iedge)=0
464  grid%CELLS(2,iedge)=0
465  lfreeedges=.true.
466  END DO
467  cycle
468  END IF
469  IF(gridman_dbg)
470  w WRITE(gridman_unit,*) " IC, IX1, IX2, IY1, IY2 ",
471  w ic,ix1,ix2,iy1,iy2
472 
473  DO iy=iy1,iy2
474 C
475 C ----2---- ----2----
476 C | | "East" "West" | |
477 C 3 1<------------------ 3 1
478 C | | edge edge | |
479 C ----4---- ----4----
480 C NXCUT(1,IC) NXCUT(2,IC)
481 C
482 C MODIFY CELLS
483  icell1=ix1-lx+1+nx*(iy-ly) !CELL ( IX1,IY )
484  icell2=ix2-lx+1+nx*(iy-ly) !CELL ( IX2,IY )
485  iedge=icell2+iy-ly+m !SIDE 3 OF ICELL2
486  ipoint1=grid%POINTS(1,iedge) !"LOWER" POINT OF THIS EDGE
487  ipoint2=grid%POINTS(2,iedge) !"UPPER" POINT OF THIS EDGE
488 C JOIN "RIGHT" EDGE WITH "LEFT" EDGE: SIDE 1 OF CELL NXCUT(1,IC)
489  grid%CELLS(1,iedge)=icell1 !MODUFY INDEX OF NEIGHBOURING
490  grid%CELLS(2,iedge)=icell2 !CELLS FOR THIS EDGE
491 C CHECK: CORRESPONDING COORDINATES OF CELL IX2 OF THE B2 GRID MUST BE SAME
492  x1=br(ix2,iy,4)
493  y1=bz(ix2,iy,4)
494  x2=br(ix2,iy,3)
495  y2=bz(ix2,iy,3)
496  IF(abs(x1-grid%X(1,ipoint1)).GT.
497  f gridman_tol*(abs(x1)+abs(grid%X(1,ipoint1))).OR.
498  f abs(y1-grid%X(2,ipoint1)).GT.
499  f gridman_tol*(abs(y1)+abs(grid%X(2,ipoint1))) .OR.
500  f abs(x2-grid%X(1,ipoint2)).GT.
501  f gridman_tol*(abs(x2)+abs(grid%X(1,ipoint2))).OR.
502  f abs(y2-grid%X(2,ipoint2)).GT.
503  f gridman_tol*(abs(y2)+abs(grid%X(2,ipoint2))) ) THEN
504  ierr=100
505  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
506  o "coordinates at the cut do not match each over"
507  WRITE(gridman_unit,*) "POINT1: X1,Y1 vs. X2,Y2 ",x1,y1,
508  w grid%X(1,ipoint1),grid%X(2,ipoint1)
509  WRITE(gridman_unit,*) "POINT2: X1,Y1 vs. X2,Y2 ",x2,y2,
510  w grid%X(1,ipoint2),grid%X(2,ipoint2)
511  WRITE(gridman_unit,*) "IX1, IX2, IY ",ix1,ix2,iy
512  GOTO 2000
513  END IF
514 C CORRECT MAPPING
515  grid%EDGEINDEX(3)%INDEXES(1,iedge-npol)=ix1 !for radial edges IE=IEDGE-NPOL
516 C CORRECT INDEX OF THE 2nd POINT OF SIDE 4 OF THE CELL NXCUT(1,IC)
517  iedge=icell1 !THIS MUST BE FIRST
518  grid%POINTS(2,iedge)=ipoint1 !("LOWER") POINT
519  END DO !DO IY=IY1,IY2
520 C IDENTIFY X-POINT
521  iedge=icell1+m+iy-ly
522  ipoint=grid%POINTS(2,iedge)
523  x1=grid%X(1,ipoint)
524  y1=grid%X(2,ipoint)
525  x2=grid%X(1,ipoint2)
526  y2=grid%X(2,ipoint2)
527  IF(abs(x1-x2).GT.gridman_tol*(abs(x1)+abs(x2)).OR.
528  f abs(y1-y2).GT.gridman_tol*(abs(x1)+abs(x2))) THEN
529  ierr=-100
530  WRITE(gridman_unit,*) "WARNING from GRIDMAN_CARRE2GRID: ",
531  o "could not identtify X-point"
532  WRITE(gridman_unit,*) "POINT1: X1,Y1 ",x1,y1
533  WRITE(gridman_unit,*) "POINT2: X2,Y2 ",x2,y2
534  END IF
535  grid%POINTS(2,iedge)=ipoint2
536 C
537 C CORRECTING POINT INDICES AT THE X-POINT FOR NXCUT(1,IC)
538 C (AT NXCUT(2,IC) THE POINT WAS NOT CHANGED)
539 C UPPER END OF THE CUT
540 C
541 C -----2------- --------
542 C |
543 C ICELL1+NX 1 X-POINT
544 C |/
545 C -----4------- --------
546 C /|
547 C ICELL1 1 | ICELL1+1 , ICELL IS THE CELL ( IX1, IX2 )
548 C / |
549 C \"JOINT" EDGE
550 C
551  icell1=icell1+nx !MOVE UP
552  iedge=icell1+m+iy-ly+1 !SIDE 1 OF THE CELL ICELL1
553  grid%POINTS(1,iedge)=ipoint2 !"LOWER" POINT IS REPLACED BY THE X-POINT
554  iedge=icell1 !SIDE 4 OF THE CELL ICELL1
555  grid%POINTS(2,iedge)=ipoint2 !"RIGHT" POINT IS REPLACED BY THE X-POINT
556  IF(ix1.LT.nxl) THEN
557  icell1=icell1+1 !MOVE RIGHT
558  iedge=icell1 !SIDE 4 OF THE CELL ICELL1
559  grid%POINTS(1,iedge)=ipoint2 !"LEFT" POINT IS REPLACED BY THE X-POINT
560  END IF
561  IF(iy1.GT.ly) THEN
562  ierr=400
563  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
564  o "NYCUT(1)>MIN(IY) is not implemented"
565  WRITE(gridman_unit,*) " NYCUT(1) ",nycut(1)
566  GOTO 2000
567  END IF
568  END DO !DO IC=1,NNCUT
569 
570  IF(lfreeedges) THEN
571 C TO CLEAN UP
572  CALL eliminate_free_edges(ierr)
573  IF(ierr.NE.0) GOTO 2000
574  END IF
575 C
576 C ASSIGN MARKERS OF THE OPEN SURFACES:
577 C -1: north boundary; -2: west; -3: east; -4: north.
578 C
579  IF(nniso.GT.0) THEN
580  CALL eliminate_isolated_cells(ierr)
581  IF(ierr.NE.0) GOTO 2000
582  END IF
583 
584 C
585 C ELIMINATE ISOLATED CELLS
586 C
587  CALL assign_surface_markers(ierr)
588  IF(ierr.NE.0) GOTO 2000
589 
590 C
591 C SHIFT CELL INDICES FOR CUTS IN THE SENSE OF EIRENE !!!!!!!
592 C
593  IF(PRESENT(leir)) THEN
594  IF(leir) THEN
595  DO ic=1,nncut
596  IF(lfort30) THEN
597  ix1=nxcut30(1,ic)+ic-1
598  ELSE
599  ix1=nxcut(ic)+ic-1
600  END IF
601  DO icell=1,grid%CELLINDEX(1)%NELEMENTS
602  ix=grid%CELLINDEX(1)%INDEXES(1,icell)
603  IF(ix.GT.ix1) grid%CELLINDEX(1)%INDEXES(1,icell)=ix+1
604  END DO
605  END DO
606  DO ic=1,nniso
607  IF(lfort30) THEN
608  ix1=nxiso30(1,ic)+ic+nncut-2
609  ELSE
610  ix1=nxiso(ic)+ic+nncut-2
611  END IF
612  DO icell=1,grid%CELLINDEX(1)%NELEMENTS
613  ix=grid%CELLINDEX(1)%INDEXES(1,icell)
614  IF(ix.GT.ix1) grid%CELLINDEX(1)%INDEXES(1,icell)=ix+1
615  END DO
616  END DO
617  END IF !IF(LEIR) THEN
618  END IF !IF(PRESENT(LEIR))
619 
620 C UNITS
621  grid%UNIT2SI=1.0_gridman_dp
622  grid%UNITS='METER'
623 C DEFAULT DESCRIPTION
624  grid%DESCRIPTION='CARRE grid (B2 plasma grid)'
625  grid%CELLINDEX(1)%DESCRIPTION=
626  = "Indexes IX, IY of the plasma grid cells"
627  grid%CELLINDEX(1)%COLUMNS(1)="CARRE_CELL_IX"
628  grid%CELLINDEX(1)%COLUMNS(2)="CARRE_CELL_IY"
629 
630  grid%EDGEINDEX(1)%DESCRIPTION= "Markers of grid boundaries:"
631  / //" -1 is South, -2 is West, -3 is East, -4 is North"
632  grid%EDGEINDEX(1)%COLUMNS(1)="CARRE_EDGE_ITAG"
633 
634  grid%EDGEINDEX(2)%DESCRIPTION=
635  = "Indexes IX, IY of the poloidal edges of plasma grid"
636  grid%EDGEINDEX(2)%COLUMNS(1)="CARRE_POL_FACE_IX"
637  grid%EDGEINDEX(2)%COLUMNS(2)="CARRE_POL_FACE_IY"
638 
639  grid%EDGEINDEX(3)%DESCRIPTION=
640  = "Indexes IX, IY of the radial edges of plasma grid"
641  grid%EDGEINDEX(3)%COLUMNS(1)="CARRE_RAD_FACE_IX"
642  grid%EDGEINDEX(3)%COLUMNS(2)="CARRE_RAD_FACE_IY"
643 
644  IF(gridman_check) THEN
645  CALL gridman_grid2d_check(grid,res,ierr)
646  IF(res.NE.0.OR.ierr.GT.0) THEN
647  ierr=400
648  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
649  w "resulting grid object is incorrect"
650  END IF
651  END IF
652 
653  IF(gridman_dbg)
654  w WRITE(gridman_unit,*) "GRIDMAN_CARRE2GRID finished"
655  RETURN
656 
657  1000 WRITE(gridman_unit,*) "GRIDMAN_CARRE2GRID terminated"
658  2000 CALL gridman_grid_deallocate(grid,ierr0)
659  RETURN
660 
661  CONTAINS
662 
663 C
664  SUBROUTINE eliminate_free_edges(IERR)
666  INTEGER,INTENT(OUT) :: IERR
667  TYPE(gridman_grid) :: GRID_TMP
668  INTEGER :: IERR0
669 
670  CALL gridman_grid_copy(grid_tmp,grid,ierr)
671  IF(ierr.NE.0) RETURN
672 
673  CALL gridman_grid_remove_free_edges(grid,grid_tmp,ierr)
674  CALL gridman_grid_deallocate(grid_tmp,ierr0)
675 
676  END SUBROUTINE eliminate_free_edges
677 
678 C
679  SUBROUTINE eliminate_isolated_cells(IERR)
681  INTEGER,INTENT(OUT) :: IERR
682  LOGICAL :: LTAKE(grid%ncells)
683  TYPE(gridman_grid) :: GRID_TMP
684  INTEGER(GRIDMAN_SP) :: ISO,ICELL,IX,IY
685  INTEGER :: IERR0
686 
687  ierr=0
688 
689  IF(nniso.LT.1) RETURN
690 
691  ltake=.true.
692 
693  IF(lfort30) THEN
694  DO iso=1,nniso
695  IF(gridman_dbg) THEN
696  WRITE(gridman_unit,*) " NXISO ",nxiso30(1,iso),nxiso30(2,iso)
697  WRITE(gridman_unit,*) " NYISO ",nyiso30(1,iso),nyiso30(2,iso)
698  END IF
699  DO icell=1,grid%NCELLS
700  ix=grid%CELLINDEX(1)%INDEXES(1,icell)
701  iy=grid%CELLINDEX(1)%INDEXES(2,icell)
702  IF(ix.GE.nxiso30(1,iso).AND.ix.LE.nxiso30(2,iso)+1.AND.
703  f iy.GE.nyiso30(1,iso).AND.iy.LE.nyiso30(2,iso)+1) THEN
704  ltake(icell)=.false.
705  END IF
706  END DO
707  END DO
708  ELSE
709  DO iso=1,nniso
710  IF(gridman_dbg) THEN
711  WRITE(gridman_unit,*) " NXISO ",nxiso(iso)
712  WRITE(gridman_unit,*) " NYISO ",nyiso(iso)
713  END IF
714  DO icell=1,grid%NCELLS
715  ix=grid%CELLINDEX(1)%INDEXES(1,icell)
716  iy=grid%CELLINDEX(1)%INDEXES(2,icell)
717  IF(ix.EQ.nxiso(iso).OR.ix.EQ.nxiso(iso)+1) THEN
718  ltake(icell)=.false.
719  END IF
720  END DO
721  END DO
722  END IF !IF(LFORT30) THEN
723 
724  CALL gridman_grid_copy(grid_tmp,grid,ierr)
725  IF(ierr.NE.0) GOTO 50
726  CALL gridman_grid_eliminate_cells(grid,grid_tmp,ltake,ierr)
727 
728  50 CALL gridman_grid_deallocate(grid_tmp,ierr0)
729 
730  END SUBROUTINE eliminate_isolated_cells
731 
732 C
733  SUBROUTINE assign_surface_markers(IERR)
734  INTEGER(GRIDMAN_SP) :: ICELL1,ICELL2,IEDGE,IE,N,IX,IY,
735  i ixy(grid%NEDGES)
736  LOGICAL :: LPOL(grid%nedges)
737  INTEGER,INTENT(OUT) :: IERR
738 
739  ierr=0
740 
741 C Count the number of elements
742  n=0
743  DO iedge=1,grid%NEDGES
744  icell1=grid%CELLS(1,iedge)
745  icell2=grid%CELLS(2,iedge)
746  IF(icell1.LT.1.AND.icell2.LT.1) cycle
747  IF(icell1.LT.1.OR.icell2.LT.1) n=n+1
748  END DO
749 
750 C Allocate index
751  CALL gridman_index_allocate(grid%EDGEINDEX(1),1,n,ierr)
752  IF(ierr.NE.0) THEN
753  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
754  w "cannot allocate index for surface markers"
755  ierr=200
756  RETURN
757  END IF
758 
759 C IY of poloidal edges
760  lpol=.false.
761  DO ie=1,grid%EDGEINDEX(2)%NELEMENTS
762  iedge=grid%EDGEINDEX(2)%INDEXES(0,ie)
763  lpol(iedge)=.true.
764  ixy(iedge)=grid%EDGEINDEX(2)%INDEXES(2,ie) !IY
765  END DO
766 C IX of radial edges
767  DO ie=1,grid%EDGEINDEX(3)%NELEMENTS
768  iedge=grid%EDGEINDEX(3)%INDEXES(0,ie)
769  IF(lpol(iedge)) THEN
770  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
771  w "internal error"
772  WRITE(gridman_unit,*) "Same edge is both poloidal and radial"
773  WRITE(gridman_unit,*) " IEDGE, IY ",iedge,ixy(iedge)
774  ierr=400
775  RETURN
776  END IF
777  ixy(iedge)=grid%EDGEINDEX(3)%INDEXES(1,ie) !IX
778  END DO
779 
780 C Assign indexes
781  ie=0
782  DO iedge=1,grid%NEDGES
783  icell1=grid%CELLS(1,iedge)
784  icell2=grid%CELLS(2,iedge)
785  IF(icell1.LT.1.OR.icell2.LT.1) THEN
786 C OPEN EDGE
787  IF(icell1.GT.0) THEN
788  ix=grid%CELLINDEX(1)%INDEXES(1,icell1)
789  iy=grid%CELLINDEX(1)%INDEXES(2,icell1)
790  ELSE
791  IF(icell2.LT.1) cycle
792  ix=grid%CELLINDEX(1)%INDEXES(1,icell2)
793  iy=grid%CELLINDEX(1)%INDEXES(2,icell2)
794  END IF
795  ie=ie+1
796  IF(ie.GT.n) THEN
797  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
798  w "internal error"
799  WRITE(gridman_unit,*)
800  w "Index is out of range in ASSIGN_SURFACE_MARKERS"
801  WRITE(gridman_unit,*) " IE, N ",ie,n
802  ierr=400
803  RETURN
804  END IF
805  grid%EDGEINDEX(1)%INDEXES(0,ie)=iedge
806  IF(lpol(iedge)) THEN
807 C POLOIDAL EDGE
808  IF(ixy(iedge).LT.iy) THEN
809 C South edge
810  grid%EDGEINDEX(1)%INDEXES(1,ie)=-1
811  ELSE
812 C North edge
813  grid%EDGEINDEX(1)%INDEXES(1,ie)=-4
814  END IF
815  ELSE
816 C RADIAL EDGE
817  IF(ixy(iedge).LT.ix) THEN
818 C West edge
819  grid%EDGEINDEX(1)%INDEXES(1,ie)=-2
820  ELSE
821 C East edge
822  grid%EDGEINDEX(1)%INDEXES(1,ie)=-3
823  END IF
824  END IF ! IF(LPOL(IEDGE))
825  END IF ! IF(ICELL1.LT.1.OR.ICELL2.LT.1)
826  END DO ! DO IEDGE=1,GRID%NEDGES
827 
828  IF(ie.NE.n) THEN
829  WRITE(gridman_unit,*) "ERROR in GRIDMAN_CARRE2GRID: ",
830  w "internal error"
831  WRITE(gridman_unit,*) "Mismatch in ASSIGN_SURFACE_MARKERS"
832  WRITE(gridman_unit,*) " IE, N ",ie,n
833  ierr=400
834  RETURN
835  END IF
836 
837 
838  END SUBROUTINE assign_surface_markers
839 
840 
841  END SUBROUTINE gridman_carre2grid
842 
843 C*********************************************************************************
844 C> Read B2 (CARRE, SONNET) grid from fort.30, return GRID object
845 C*********************************************************************************
846 C>
847 C> GRIDMAN_CARRE_READ30_GRID is a shell for GRIDMAN_CARRE_READ30_ARRAY
848 C>
849 C> GRIDMAN_CARRE_READ30_ARRAY and GRIDMAN_CARRE_READ30_GRID
850 C> are combined in the interface GRIDMAN_CARRE_READ30
851 C>
852  SUBROUTINE gridman_carre_read30_grid(GRID,FNAME,IERR,LEIR)
854  u gridman_unit,gridman_dbg,gridman_length
856  u gridman_carre2grid_fort30,
858  IMPLICIT NONE
859 C> Resulting grid object
860  TYPE(gridman_grid) :: GRID
861 C> String containing the name of input file
862  CHARACTER(*),INTENT(IN) :: FNAME
863 C> Error code
864  INTEGER,INTENT(OUT) :: IERR
865 C> Switch to produce cell indexes (of the original plasma grid)
866 C> in the sense of EIRENE
867 C>
868 C> If .TRUE. then inices in the index mapping of cells are shifted
869 C> in the presence if cuts in the way it is done in EIRENE
870  LOGICAL,OPTIONAL :: LEIR
871 
872  CHARACTER(GRIDMAN_LENGTH) :: SONNETFILE
873  INTEGER(GRIDMAN_SP) :: NX,NY,NNCUT,NNISO
874  INTEGER :: IERR0
875  INTEGER(GRIDMAN_SP),ALLOCATABLE :: NXCUT(:,:),NYCUT(:,:),
876  i nxiso(:,:),nyiso(:,:)
877  REAL(GRIDMAN_DP),ALLOCATABLE :: BR(:,:,:),BZ(:,:,:)
878 
879  IF(gridman_dbg)
880  w WRITE(gridman_unit,*) "Starting GRIDMAN_CARRE_READ30_GRID"
881 
882  ierr=0
883 
884  CALL gridman_carre_read30_array(fname,sonnetfile,nx,ny,
885  s nncut,nxcut,nycut,
886  s nniso,nxiso,nyiso,br,bz,ierr)
887  IF(ierr.NE.0) THEN
888  WRITE(gridman_unit,*) "GRIDMAN_CARRE_READ30_GRID terminated"
889  RETURN
890  END IF
891 
892  CALL gridman_carre2grid_fort30(grid,1_gridman_sp,1_gridman_sp,
893  s nx,ny,nncut,nxcut,nycut,
894  s nniso,nxiso,nyiso,br,bz,
895  s 1_gridman_sp,1_gridman_sp,
896  s nx,ny,ierr,leir)
897 
898  IF(ALLOCATED(nxcut)) DEALLOCATE(nxcut)
899  IF(ALLOCATED(nycut)) DEALLOCATE(nycut)
900  IF(ALLOCATED(nxiso)) DEALLOCATE(nxiso)
901  IF(ALLOCATED(nyiso)) DEALLOCATE(nyiso)
902  IF(ALLOCATED(br)) DEALLOCATE(br)
903  IF(ALLOCATED(bz)) DEALLOCATE(bz)
904 
905  IF(ierr.NE.0) THEN
906  CALL gridman_grid_deallocate(grid,ierr0)
907  WRITE(gridman_unit,*) "GRIDMAN_CARRE_READ30_GRID terminated"
908  RETURN
909  END IF
910 
911  grid%DESCRIPTION=sonnetfile
912 
913  IF(gridman_dbg)
914  w WRITE(gridman_unit,*) "GRIDMAN_CARRE_READ30_GRID finished"
915 
916  END SUBROUTINE gridman_carre_read30_grid
917 
918 
919 C*********************************************************************************
920 C> Read CARRE grid in SONNET format, return GRID object
921 C*********************************************************************************
922 C>
923 C> GRIDMAN_CARRE_READSONNET_ARRAY and GRIDMAN_CARRE_READSONNET_GRID
924 C> are combined in the interface GRIDMAN_CARRE_READSONNET
925 C>
926  SUBROUTINE gridman_carre_readsonnet_grid(GRID,FNAME,IERR,LEIR)
928  u gridman_unit,gridman_dbg
931  IMPLICIT NONE
932 C> Resulting grid object
933  TYPE(gridman_grid) :: GRID
934 C> String containing the name of input file
935  CHARACTER(*),INTENT(IN) :: FNAME
936 C> Error code
937  INTEGER,INTENT(OUT) :: IERR
938 C> Switch to produce cell index mapping in the sense of EIRENE
939 C>
940 C> If .TRUE. then inices in the index mapping of cells are shifted
941 C> in the presence if cuts in the way it is done in EIRENE
942  LOGICAL,OPTIONAL :: LEIR
943 
944  INTEGER(GRIDMAN_SP) :: NX,NY,NNCUT,NNISO
945  INTEGER :: IERR0
946  INTEGER(GRIDMAN_SP),ALLOCATABLE :: NXCUT(:),NYCUT(:),
947  i nxiso(:),nyiso(:)
948  REAL(GRIDMAN_DP),ALLOCATABLE :: BR(:,:,:),BZ(:,:,:),
949  r bc(:,:,:),pit(:,:)
950  REAL(GRIDMAN_DP) :: RBT
951 
952  IF(gridman_dbg)
953  w WRITE(gridman_unit,*) "Starting GRIDMAN_CARRE_READSONNET_GRID"
954 
955  ierr=0
956 
957  CALL gridman_carre_readsonnet_array(fname,nx,ny,
958  s nncut,nxcut,nycut,
959  s nniso,nxiso,nyiso,
960  s rbt,br,bz,pit,bc,ierr)
961  IF(ierr.NE.0) THEN
962  WRITE(gridman_unit,*) "GRIDMAN_CARRE_READSONNET_GRID terminated"
963  RETURN
964  END IF
965 
966  CALL gridman_carre2grid(grid,0_gridman_sp,0_gridman_sp,
967  s nx+1,ny+1,nncut,
968  s nxcut,nycut,nniso,nxiso,nyiso,
969  s br,bz,1_gridman_sp,1_gridman_sp,
970  s nx,ny,ierr,leir)
971 
972  IF(ALLOCATED(nxcut)) DEALLOCATE(nxcut)
973  IF(ALLOCATED(nycut)) DEALLOCATE(nycut)
974  IF(ALLOCATED(nxiso)) DEALLOCATE(nxiso)
975  IF(ALLOCATED(nyiso)) DEALLOCATE(nyiso)
976  IF(ALLOCATED(br)) DEALLOCATE(br)
977  IF(ALLOCATED(bz)) DEALLOCATE(bz)
978  IF(ALLOCATED(bc)) DEALLOCATE(bc)
979  IF(ALLOCATED(pit)) DEALLOCATE(pit)
980 
981  IF(ierr.NE.0) THEN
982  CALL gridman_grid_deallocate(grid,ierr0)
983  WRITE(gridman_unit,*) "GRIDMAN_CARRE_READSONNET_GRID terminated"
984  RETURN
985  END IF
986 
987  IF(gridman_dbg)
988  w WRITE(gridman_unit,*) "GRIDMAN_CARRE_READSONNET_GRID finished"
989 
990  END SUBROUTINE gridman_carre_readsonnet_grid
subroutine gridman_carre2grid(GRID, LXB, LYB, NXB, NYB, NNCUT, NXCUT, NYCUT, NNISO, NXISO, NYISO, BR, BZ, LX, LY, NXL, NYL, IERR, LEIR)
Convert B2 (CARRE, SONNET) grid into GRIDMAN_GRID grid object.
Definition: carre.f:72
real(gridman_dp), save, public gridman_tol
Tolerance parameter which is used to compare two real numbers.
Definition: gridman.f:127
subroutine gridman_grid_eliminate_cells(GRID_NEW, GRID, LTAKE, IERR)
Eliminate cells from GRIDMAN_GRID object.
Definition: grid2.f:544
integer, parameter, public gridman_dp
Kind parameter for real numbers.
Definition: gridman.f:93
Explicit interfaces to GRIDMAN subroutines and functions.
Definition: gridman.f:251
subroutine gridman_grid_allocate(GRID, TYPE, NEDGES, NPOINTS, NCELLS, IERR, NEDGEINDEX, NCELLINDEX)
Allocate GRIDMAN_GRID object.
Definition: grid1.f:30
subroutine gridman_carre_read30_array(FNAME, STITLE, NX, NY, NNCUT, NXCUT, NYCUT, NNISO, NXISO, NYISO, BR, BZ, IERR)
Read B2 (CARRE, SONNET) grid from fort.30, return separate arrays.
Definition: carre.f:54
subroutine gridman_grid2d_check(GRID, RES, IERR)
Check correctness of the 2D grid object.
Definition: grid2d.f:37
Data-type which describes a grid as a set of edges, methods in grid.f.
Definition: gridman.f:168
subroutine gridman_grid_deallocate(GRID, IERR)
Deallocate grid object.
Definition: grid1.f:184
subroutine gridman_index_allocate(INDEX, NINDEX, NELEMENTS, IERR)
Allocate index object.
Definition: index.f:28
subroutine gridman_grid_copy(GRID2, GRID1, IERR)
Create a copy of the grid object.
Definition: grid1.f:981
subroutine gridman_carre_readsonnet_array(FNAME, NX, NY, NNCUT, NXCUT, NYCUT, NNISO, NXISO, NYISO, RBT, BR, BZ, PIT, BC, IERR)
Read carre grid in SONNET format, return separate arrays.
Definition: carre.f:282
subroutine gridman_carre_readsonnet_grid(GRID, FNAME, IERR, LEIR)
Read CARRE grid in SONNET format, return GRID object.
Definition: carre.f:927
Definition of data types, global constants and variables.
Definition: gridman.f:83
subroutine gridman_carre_read30_grid(GRID, FNAME, IERR, LEIR)
Read B2 (CARRE, SONNET) grid from fort.30, return GRID object.
Definition: carre.f:853
subroutine gridman_grid_remove_free_edges(GRID_NEW, GRID, IERR)
Remove edges which do not belong to any cell from the GRIDMAN_GRID object.
Definition: grid2.f:999
integer, parameter, public gridman_sp
Kind parameter for integer numbers.
Definition: gridman.f:95