ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
coarseUtils.F90
Go to the documentation of this file.
1 module coarseutils
2 
3 contains
4 
5  subroutine createcoarseblocks(level)
6  !
7  ! createCoarseBlocks creates the block data structure for the
8  ! given coarse grid from the 1 level finer grid. Only direct
9  ! info is created, like owned coordinates, block sizes and
10  ! subface info. Indirect info, like face normals, volumes, wall
11  ! distances, etc. Are created later on. That info can be created
12  ! independent of the finer grid.
13  !
14  use constants
15  use block
18  use coarseninginfo
19  use utils, only: terminate
20  implicit none
21  !
22  ! Subroutine arguments.
23  !
24  integer(kind=intType), intent(in) :: level
25  !
26  ! Local variables.
27  !
28  integer :: ierr
29 
30  integer(kind=intType) :: i, j, k, ii, jj, kk, n1to1
31  integer(kind=intType) :: nn, mm, iil, jjl, kkl, il, jl, kl
32  integer(kind=intType) :: iBeg, jBeg, kBeg, iEnd, jEnd, kEnd
33  integer(kind=intType) :: l1, L2, l3
34  integer(kind=intType) :: levm1, donorOffset, fact
35 
36  integer(kind=intType), dimension(:), allocatable :: imap, iimap
37  integer(kind=intType), dimension(:), allocatable :: jmap, jjmap
38  integer(kind=intType), dimension(:), allocatable :: kmap, kkmap
39 
40  integer(kind=intType), dimension(:), pointer :: dfine
41 
42  logical, dimension(:), pointer :: iCo, jCo, kCo
43 
44  ! Store the finer grid level in levm1
45 
46  levm1 = level - 1
47 
48  ! Determine the total number of 1 to 1 block faces on this
49  ! processor for the fine grid level.
50 
51  nsubface1to1 = 0
52  do nn = 1, ndom
53  nsubface1to1 = nsubface1to1 + flowdoms(nn, levm1, 1)%n1to1
54  end do
55 
56  ! Allocate the memory for subface1to1 and coarseInfo.
57 
58  allocate (subface1to1(nsubface1to1), coarseinfo(ndom), stat=ierr)
59  if (ierr /= 0) &
60  call terminate("createCoarseBlocks", &
61  "Memory allocation failure for subface1to1 &
62  &and coarseInfo")
63 
64  ! Loop over the number of blocks on this processor.
65 
66  n1to1 = 0
67  domains: do nn = 1, ndom
68 
69  ! If levm1 is 1, i.e. the finest grid, set the coarsenings to
70  ! regular, as this is a regular grid.
71 
72  if (levm1 == 1) then
73  flowdoms(nn, levm1, 1)%iCoarsened = regular
74  flowdoms(nn, levm1, 1)%jCoarsened = regular
75  flowdoms(nn, levm1, 1)%kCoarsened = regular
76  end if
77 
78  ! Copy cgns block id, blockIsMoving and addGridVelocities.
79  ! Although identical to the values on the finest grid level,
80  ! it is copied anyway for consistency reasons.
81 
82  flowdoms(nn, level, 1)%cgnsBlockID = &
83  flowdoms(nn, levm1, 1)%cgnsBlockID
84  flowdoms(nn, level, 1)%blockIsMoving = &
85  flowdoms(nn, levm1, 1)%blockIsMoving
86  flowdoms(nn, level, 1)%addGridVelocities = &
87  flowdoms(nn, levm1, 1)%addGridVelocities
88 
89  ! Store the number of fine nodes a bit easier and allocate the
90  ! memory for the logicals iCo, jCo and kCo and for coarseIs1to1.
91 
92  iil = flowdoms(nn, levm1, 1)%il
93  jjl = flowdoms(nn, levm1, 1)%jl
94  kkl = flowdoms(nn, levm1, 1)%kl
95 
96  ii = flowdoms(nn, levm1, 1)%n1to1
97  allocate (flowdoms(nn, levm1, 1)%iCo(iil), &
98  flowdoms(nn, levm1, 1)%jCo(jjl), &
99  flowdoms(nn, levm1, 1)%kCo(kkl), &
100  coarseinfo(nn)%coarseIs1to1(ii), stat=ierr)
101  if (ierr /= 0) &
102  call terminate("createCoarseBlocks", &
103  "Memory allocation failure for iCo, jCo, kCo &
104  &and coarseIs1to1")
105 
106  ! Initialize iCo, jCo and kCo such that the block boundaries
107  ! remain, but all other nodes can disappear. Set pointers to
108  ! make the code more readable.
109 
110  ico => flowdoms(nn, levm1, 1)%iCo
111  jco => flowdoms(nn, levm1, 1)%jCo
112  kco => flowdoms(nn, levm1, 1)%kCo
113 
114  ico = .false.; ico(1) = .true.; ico(iil) = .true.
115  jco = .false.; jco(1) = .true.; jco(jjl) = .true.
116  kco = .false.; kco(1) = .true.; kco(kkl) = .true.
117 
118  ! Loop over the subfaces to keep their boundaries. Also internal
119  ! block boundaries are kept.
120 
121  do i = 1, flowdoms(nn, levm1, 1)%nSubface
122  ico(flowdoms(nn, levm1, 1)%inBeg(i)) = .true.
123  ico(flowdoms(nn, levm1, 1)%inEnd(i)) = .true.
124 
125  jco(flowdoms(nn, levm1, 1)%jnBeg(i)) = .true.
126  jco(flowdoms(nn, levm1, 1)%jnEnd(i)) = .true.
127 
128  kco(flowdoms(nn, levm1, 1)%knBeg(i)) = .true.
129  kco(flowdoms(nn, levm1, 1)%knEnd(i)) = .true.
130  end do
131 
132  ! Create the coarser grid in i-direction. In case the fine grid
133  ! was already a nonregular coarsened block, start at the
134  ! opposite boundary. For a regular grid simply start at the
135  ! left boundary.
136 
137  if (flowdoms(nn, levm1, 1)%iCoarsened == leftstarted) then
138 
139  flowdoms(nn, level, 1)%iCoarsened = rightstarted
140 
141  do i = (iil - 1), 2, -1
142  if (.not. ico(i + 1)) ico(i) = .true.
143  end do
144 
145  else
146 
147  flowdoms(nn, level, 1)%iCoarsened = leftstarted
148 
149  do i = 2, (iil - 1)
150  if (.not. ico(i - 1)) ico(i) = .true.
151  end do
152 
153  end if
154 
155  ! Create the coarser grid in j-direction. Same story as in
156  ! i-direction.
157 
158  if (flowdoms(nn, levm1, 1)%jCoarsened == leftstarted) then
159 
160  flowdoms(nn, level, 1)%jCoarsened = rightstarted
161 
162  do j = (jjl - 1), 2, -1
163  if (.not. jco(j + 1)) jco(j) = .true.
164  end do
165 
166  else
167 
168  flowdoms(nn, level, 1)%jCoarsened = leftstarted
169 
170  do j = 2, (jjl - 1)
171  if (.not. jco(j - 1)) jco(j) = .true.
172  end do
173 
174  end if
175 
176  ! Create the coarser grid in k-direction. Same story as in
177  ! i- and j-direction.
178 
179  if (flowdoms(nn, levm1, 1)%kCoarsened == leftstarted) then
180 
181  flowdoms(nn, level, 1)%kCoarsened = rightstarted
182 
183  do k = (kkl - 1), 2, -1
184  if (.not. kco(k + 1)) kco(k) = .true.
185  end do
186 
187  else
188 
189  flowdoms(nn, level, 1)%kCoarsened = leftstarted
190 
191  do k = 2, (kkl - 1)
192  if (.not. kco(k - 1)) kco(k) = .true.
193  end do
194 
195  end if
196 
197  ! Determine the number of points in each direction for the
198  ! coarse grid.
199 
200  il = 0
201  do i = 1, iil
202  if (ico(i)) il = il + 1
203  end do
204 
205  jl = 0
206  do j = 1, jjl
207  if (jco(j)) jl = jl + 1
208  end do
209 
210  kl = 0
211  do k = 1, kkl
212  if (kco(k)) kl = kl + 1
213  end do
214 
215  ! Store the number of nodes and cells in the three directions
216  ! and the dimensions for the halo based quantities.
217 
218  flowdoms(nn, level, 1)%il = il
219  flowdoms(nn, level, 1)%jl = jl
220  flowdoms(nn, level, 1)%kl = kl
221 
222  flowdoms(nn, level, 1)%nx = il - 1
223  flowdoms(nn, level, 1)%ny = jl - 1
224  flowdoms(nn, level, 1)%nz = kl - 1
225 
226  flowdoms(nn, level, 1)%ie = il + 1
227  flowdoms(nn, level, 1)%je = jl + 1
228  flowdoms(nn, level, 1)%ke = kl + 1
229 
230  flowdoms(nn, level, 1)%ib = il + 2
231  flowdoms(nn, level, 1)%jb = jl + 2
232  flowdoms(nn, level, 1)%kb = kl + 2
233 
234  ! If the coarsening was regular, i.e. if the number of fine
235  ! grid cells is twice the number of coarse grid cells, reset
236  ! the coarsening to regular.
237 
238  if (flowdoms(nn, levm1, 1)%nx == 2 * flowdoms(nn, level, 1)%nx) &
239  flowdoms(nn, level, 1)%iCoarsened = regular
240  if (flowdoms(nn, levm1, 1)%ny == 2 * flowdoms(nn, level, 1)%ny) &
241  flowdoms(nn, level, 1)%jCoarsened = regular
242  if (flowdoms(nn, levm1, 1)%nz == 2 * flowdoms(nn, level, 1)%nz) &
243  flowdoms(nn, level, 1)%kCoarsened = regular
244  !
245  ! The variables, which control the restriction to and the
246  ! interpolation from the coarser grid level.
247  !
248  ! Allocate the memory.
249 
250  i = flowdoms(nn, level, 1)%ie
251  j = flowdoms(nn, level, 1)%je
252  k = flowdoms(nn, level, 1)%ke
253 
254  allocate (flowdoms(nn, level, 1)%mgIFine(1:i, 2), &
255  flowdoms(nn, level, 1)%mgJFine(1:j, 2), &
256  flowdoms(nn, level, 1)%mgKFine(1:k, 2), &
257  flowdoms(nn, level, 1)%mgIWeight(2:il), &
258  flowdoms(nn, level, 1)%mgJWeight(2:jl), &
259  flowdoms(nn, level, 1)%mgKWeight(2:kl), &
260  flowdoms(nn, levm1, 1)%mgICoarse(2:iil, 2), &
261  flowdoms(nn, levm1, 1)%mgJCoarse(2:jjl, 2), &
262  flowdoms(nn, levm1, 1)%mgKCoarse(2:kkl, 2), stat=ierr)
263  if (ierr /= 0) &
264  call terminate("createCoarseBlocks", &
265  "Memory allocation failure for interpolation &
266  &variables")
267 
268  ! Set the halo values to the halo indices of the fine level.
269 
270  flowdoms(nn, level, 1)%mgIFine(1, 1:2) = (/0, 1/)
271  flowdoms(nn, level, 1)%mgJFine(1, 1:2) = (/0, 1/)
272  flowdoms(nn, level, 1)%mgKFine(1, 1:2) = (/0, 1/)
273 
274  flowdoms(nn, level, 1)%mgIFine(i, 1) = flowdoms(nn, levm1, 1)%ie
275  flowdoms(nn, level, 1)%mgIFine(i, 2) = flowdoms(nn, levm1, 1)%ib
276  flowdoms(nn, level, 1)%mgJFine(j, 1) = flowdoms(nn, levm1, 1)%je
277  flowdoms(nn, level, 1)%mgJFine(j, 2) = flowdoms(nn, levm1, 1)%jb
278  flowdoms(nn, level, 1)%mgKFine(k, 1) = flowdoms(nn, levm1, 1)%ke
279  flowdoms(nn, level, 1)%mgKFine(k, 2) = flowdoms(nn, levm1, 1)%kb
280 
281  ! Determine the restriction variables in i-direction.
282 
283  ii = 2
284  do i = 2, iil
285  if (ico(i)) then
286  if (ico(i - 1)) then
287  flowdoms(nn, level, 1)%mgIFine(ii, 1) = i
288  flowdoms(nn, level, 1)%mgIFine(ii, 2) = i
289  flowdoms(nn, level, 1)%mgIWeight(ii) = half
290  else
291  flowdoms(nn, level, 1)%mgIFine(ii, 1) = i - 1
292  flowdoms(nn, level, 1)%mgIFine(ii, 2) = i
293  flowdoms(nn, level, 1)%mgIWeight(ii) = one
294  end if
295  ii = ii + 1
296  end if
297  end do
298 
299  ! Determine the restriction variables in j-direction.
300 
301  jj = 2
302  do j = 2, jjl
303  if (jco(j)) then
304  if (jco(j - 1)) then
305  flowdoms(nn, level, 1)%mgJFine(jj, 1) = j
306  flowdoms(nn, level, 1)%mgJFine(jj, 2) = j
307  flowdoms(nn, level, 1)%mgJWeight(jj) = half
308  else
309  flowdoms(nn, level, 1)%mgJFine(jj, 1) = j - 1
310  flowdoms(nn, level, 1)%mgJFine(jj, 2) = j
311  flowdoms(nn, level, 1)%mgJWeight(jj) = one
312  end if
313  jj = jj + 1
314  end if
315  end do
316 
317  ! Determine the restriction variables in k-direction.
318 
319  kk = 2
320  do k = 2, kkl
321  if (kco(k)) then
322  if (kco(k - 1)) then
323  flowdoms(nn, level, 1)%mgKFine(kk, 1) = k
324  flowdoms(nn, level, 1)%mgKFine(kk, 2) = k
325  flowdoms(nn, level, 1)%mgKWeight(kk) = half
326  else
327  flowdoms(nn, level, 1)%mgKFine(kk, 1) = k - 1
328  flowdoms(nn, level, 1)%mgKFine(kk, 2) = k
329  flowdoms(nn, level, 1)%mgKWeight(kk) = one
330  end if
331  kk = kk + 1
332  end if
333  end do
334 
335  ! Determine the interpolation variables in i-direction.
336 
337  ii = 2
338  do i = 2, iil
339  if (ico(i)) then
340  if (ico(i - 1)) then
341  flowdoms(nn, levm1, 1)%mgICoarse(i, 1) = ii
342  flowdoms(nn, levm1, 1)%mgICoarse(i, 2) = ii
343  else
344  flowdoms(nn, levm1, 1)%mgICoarse(i, 1) = ii
345  flowdoms(nn, levm1, 1)%mgICoarse(i, 2) = ii + 1
346  end if
347  ii = ii + 1
348  else
349  flowdoms(nn, levm1, 1)%mgICoarse(i, 1) = ii
350  flowdoms(nn, levm1, 1)%mgICoarse(i, 2) = ii - 1
351  end if
352  end do
353 
354  ! Determine the interpolation variables in j-direction.
355 
356  jj = 2
357  do j = 2, jjl
358  if (jco(j)) then
359  if (jco(j - 1)) then
360  flowdoms(nn, levm1, 1)%mgJCoarse(j, 1) = jj
361  flowdoms(nn, levm1, 1)%mgJCoarse(j, 2) = jj
362  else
363  flowdoms(nn, levm1, 1)%mgJCoarse(j, 1) = jj
364  flowdoms(nn, levm1, 1)%mgJCoarse(j, 2) = jj + 1
365  end if
366  jj = jj + 1
367  else
368  flowdoms(nn, levm1, 1)%mgJCoarse(j, 1) = jj
369  flowdoms(nn, levm1, 1)%mgJCoarse(j, 2) = jj - 1
370  end if
371  end do
372 
373  ! Determine the interpolation variables in k-direction.
374 
375  kk = 2
376  do k = 2, kkl
377  if (kco(k)) then
378  if (kco(k - 1)) then
379  flowdoms(nn, levm1, 1)%mgKCoarse(k, 1) = kk
380  flowdoms(nn, levm1, 1)%mgKCoarse(k, 2) = kk
381  else
382  flowdoms(nn, levm1, 1)%mgKCoarse(k, 1) = kk
383  flowdoms(nn, levm1, 1)%mgKCoarse(k, 2) = kk + 1
384  end if
385  kk = kk + 1
386  else
387  flowdoms(nn, levm1, 1)%mgKCoarse(k, 1) = kk
388  flowdoms(nn, levm1, 1)%mgKCoarse(k, 2) = kk - 1
389  end if
390  end do
391  !
392  ! The coordinate mapping from fine to coarse and coarse to
393  ! fine. These are needed to determine the coarse grid subface
394  ! info.
395  !
396  ! Allocate the memory.
397 
398  allocate (imap(iil), jmap(jjl), kmap(kkl), &
399  iimap(il), jjmap(jl), kkmap(kl), stat=ierr)
400  if (ierr /= 0) &
401  call terminate("createCoarseBlocks", &
402  "Memory allocation failure for imap, etc")
403 
404  ! Set the values of imap and iimap.
405 
406  ii = 1
407  do i = 1, iil
408  if (ico(i)) then
409  imap(i) = ii
410  iimap(ii) = i
411  ii = ii + 1
412  end if
413  end do
414 
415  ! Set the values of jmap and jjmap.
416 
417  jj = 1
418  do j = 1, jjl
419  if (jco(j)) then
420  jmap(j) = jj
421  jjmap(jj) = j
422  jj = jj + 1
423  end if
424  end do
425 
426  ! Set the values of kmap and kkmap.
427 
428  kk = 1
429  do k = 1, kkl
430  if (kco(k)) then
431  kmap(k) = kk
432  kkmap(kk) = k
433  kk = kk + 1
434  end if
435  end do
436  !
437  ! The subface info. Except for the subface range all other
438  ! data can be copied. The range must be adapted and the donor
439  ! range is created later, because the coarsening info of the
440  ! donor block must be known.
441  !
442  flowdoms(nn, level, 1)%nSubface = flowdoms(nn, levm1, 1)%nSubface
443  flowdoms(nn, level, 1)%n1to1 = flowdoms(nn, levm1, 1)%n1to1
444  flowdoms(nn, level, 1)%nBocos = flowdoms(nn, levm1, 1)%nBocos
445  flowdoms(nn, level, 1)%nViscBocos = flowdoms(nn, levm1, 1)%nViscBocos
446 
447  ! Allocate the memory.
448 
449  mm = flowdoms(nn, level, 1)%nSubface
450  allocate (flowdoms(nn, level, 1)%BCType(mm), &
451  flowdoms(nn, level, 1)%BCFaceID(mm), &
452  flowdoms(nn, level, 1)%cgnsSubface(mm), &
453  flowdoms(nn, level, 1)%neighBlock(mm), &
454  flowdoms(nn, level, 1)%neighProc(mm), &
455  flowdoms(nn, level, 1)%groupNum(mm), &
456  flowdoms(nn, level, 1)%inBeg(mm), &
457  flowdoms(nn, level, 1)%jnBeg(mm), &
458  flowdoms(nn, level, 1)%knBeg(mm), &
459  flowdoms(nn, level, 1)%inEnd(mm), &
460  flowdoms(nn, level, 1)%jnEnd(mm), &
461  flowdoms(nn, level, 1)%knEnd(mm), &
462  flowdoms(nn, level, 1)%dinBeg(mm), &
463  flowdoms(nn, level, 1)%djnBeg(mm), &
464  flowdoms(nn, level, 1)%dknBeg(mm), &
465  flowdoms(nn, level, 1)%dinEnd(mm), &
466  flowdoms(nn, level, 1)%djnEnd(mm), &
467  flowdoms(nn, level, 1)%dknEnd(mm), &
468  flowdoms(nn, level, 1)%l1(mm), &
469  flowdoms(nn, level, 1)%l2(mm), &
470  flowdoms(nn, level, 1)%l3(mm), &
471  stat=ierr)
472  if (ierr /= 0) &
473  call terminate("createCoarseBlocks", &
474  "Memory allocation failure for subface info")
475 
476  ! Loop over the subfaces.
477 
478  subfaces: do mm = 1, flowdoms(nn, level, 1)%nSubface
479 
480  ! Determine the range of the coarse subface.
481 
482  flowdoms(nn, level, 1)%inBeg(mm) = &
483  imap(flowdoms(nn, levm1, 1)%inBeg(mm))
484  flowdoms(nn, level, 1)%jnBeg(mm) = &
485  jmap(flowdoms(nn, levm1, 1)%jnBeg(mm))
486  flowdoms(nn, level, 1)%knBeg(mm) = &
487  kmap(flowdoms(nn, levm1, 1)%knBeg(mm))
488 
489  flowdoms(nn, level, 1)%inEnd(mm) = &
490  imap(flowdoms(nn, levm1, 1)%inEnd(mm))
491  flowdoms(nn, level, 1)%jnEnd(mm) = &
492  jmap(flowdoms(nn, levm1, 1)%jnEnd(mm))
493  flowdoms(nn, level, 1)%knEnd(mm) = &
494  kmap(flowdoms(nn, levm1, 1)%knEnd(mm))
495 
496  ! Copy the rest of the subface info.
497 
498  flowdoms(nn, level, 1)%BCType(mm) = &
499  flowdoms(nn, levm1, 1)%BCType(mm)
500  flowdoms(nn, level, 1)%BCFaceID(mm) = &
501  flowdoms(nn, levm1, 1)%BCFaceID(mm)
502  flowdoms(nn, level, 1)%neighBlock(mm) = &
503  flowdoms(nn, levm1, 1)%neighBlock(mm)
504  flowdoms(nn, level, 1)%neighProc(mm) = &
505  flowdoms(nn, levm1, 1)%neighProc(mm)
506  flowdoms(nn, level, 1)%groupNum(mm) = &
507  flowdoms(nn, levm1, 1)%groupNum(mm)
508  flowdoms(nn, level, 1)%cgnsSubface(mm) = &
509  flowdoms(nn, levm1, 1)%cgnsSubface(mm)
510 
511  flowdoms(nn, level, 1)%l1(mm) = flowdoms(nn, levm1, 1)%l1(mm)
512  flowdoms(nn, level, 1)%l2(mm) = flowdoms(nn, levm1, 1)%l2(mm)
513  flowdoms(nn, level, 1)%l3(mm) = flowdoms(nn, levm1, 1)%l3(mm)
514 
515  ! Create some info if this is a 1 to 1 subface. This is stored
516  ! in the array subface1to1 and is needed to determine the donor
517  ! info and to check if this is still a 1 to 1 subface on the
518  ! coarse grid.
519 
520  subface_1to1: if (mm > flowdoms(nn, level, 1)%nBocos .and. &
521  mm <= (flowdoms(nn, level, 1)%nBocos &
522  + flowdoms(nn, level, 1)%n1to1)) then
523 
524  ! Update the counter.
525 
526  n1to1 = n1to1 + 1
527 
528  ! Store the range of the 1 to 1 subface.
529 
530  subface1to1(n1to1)%iBeg = flowdoms(nn, level, 1)%inBeg(mm)
531  subface1to1(n1to1)%jBeg = flowdoms(nn, level, 1)%jnBeg(mm)
532  subface1to1(n1to1)%kBeg = flowdoms(nn, level, 1)%knBeg(mm)
533  subface1to1(n1to1)%iEnd = flowdoms(nn, level, 1)%inEnd(mm)
534  subface1to1(n1to1)%jEnd = flowdoms(nn, level, 1)%jnEnd(mm)
535  subface1to1(n1to1)%kEnd = flowdoms(nn, level, 1)%knEnd(mm)
536 
537  ! Store the new range of this subface a bit easier.
538  ! Make sure that i1, etc contains the lowest index.
539 
540  ibeg = min(subface1to1(n1to1)%iBeg, subface1to1(n1to1)%iEnd)
541  jbeg = min(subface1to1(n1to1)%jBeg, subface1to1(n1to1)%jEnd)
542  kbeg = min(subface1to1(n1to1)%kBeg, subface1to1(n1to1)%kEnd)
543 
544  iend = max(subface1to1(n1to1)%iBeg, subface1to1(n1to1)%iEnd)
545  jend = max(subface1to1(n1to1)%jBeg, subface1to1(n1to1)%jEnd)
546  kend = max(subface1to1(n1to1)%kBeg, subface1to1(n1to1)%kEnd)
547 
548  ! And copy it back into subface1to1.
549 
550  subface1to1(n1to1)%iBeg = ibeg
551  subface1to1(n1to1)%jBeg = jbeg
552  subface1to1(n1to1)%kBeg = kbeg
553  subface1to1(n1to1)%iEnd = iend
554  subface1to1(n1to1)%jEnd = jend
555  subface1to1(n1to1)%kEnd = kend
556 
557  ! Store the processor and block ID of the donor.
558 
559  subface1to1(n1to1)%neighProc = &
560  flowdoms(nn, level, 1)%neighProc(mm)
561  subface1to1(n1to1)%neighBlock = &
562  flowdoms(nn, level, 1)%neighBlock(mm)
563 
564  ! Store the shorthand for the transformation matrix a
565  ! bit easier.
566 
567  l1 = flowdoms(nn, level, 1)%l1(mm)
568  l2 = flowdoms(nn, level, 1)%l2(mm)
569  l3 = flowdoms(nn, level, 1)%l3(mm)
570 
571  ! Determine the fine grid donor indices of the subface.
572  ! i-direction.
573 
574  fact = 1
575  if (l1 < 0) fact = -1
576  ii = iend - ibeg + 1
577 
578  l1 = abs(l1)
579  select case (l1)
580  case (1_inttype)
581  allocate (subface1to1(n1to1)%idfine(ii), stat=ierr)
582  dfine => subface1to1(n1to1)%idfine
583  subface1to1(n1to1)%ndi = ii
584  donoroffset = flowdoms(nn, levm1, 1)%dinBeg(mm)
585  case (2_inttype)
586  allocate (subface1to1(n1to1)%jdfine(ii), stat=ierr)
587  dfine => subface1to1(n1to1)%jdfine
588  subface1to1(n1to1)%ndj = ii
589  donoroffset = flowdoms(nn, levm1, 1)%djnBeg(mm)
590  case (3_inttype)
591  allocate (subface1to1(n1to1)%kdfine(ii), stat=ierr)
592  dfine => subface1to1(n1to1)%kdfine
593  subface1to1(n1to1)%ndk = ii
594  donoroffset = flowdoms(nn, levm1, 1)%dknBeg(mm)
595  end select
596 
597  if (ierr /= 0) &
598  call terminate("createCoarseBlocks", &
599  "Memory allocation failure for idfine")
600 
601  ii = 1
602  do i = ibeg, iend
603  jj = fact * (iimap(i) - flowdoms(nn, levm1, 1)%inBeg(mm))
604  dfine(ii) = jj + donoroffset
605  ii = ii + 1
606  end do
607 
608  ! Determine the fine grid donor indices of the subface.
609  ! j-direction.
610 
611  fact = 1
612  if (l2 < 0) fact = -1
613  ii = jend - jbeg + 1
614 
615  l2 = abs(l2)
616  select case (l2)
617  case (1_inttype)
618  allocate (subface1to1(n1to1)%idfine(ii), stat=ierr)
619  dfine => subface1to1(n1to1)%idfine
620  subface1to1(n1to1)%ndi = ii
621  donoroffset = flowdoms(nn, levm1, 1)%dinBeg(mm)
622  case (2_inttype)
623  allocate (subface1to1(n1to1)%jdfine(ii), stat=ierr)
624  dfine => subface1to1(n1to1)%jdfine
625  subface1to1(n1to1)%ndj = ii
626  donoroffset = flowdoms(nn, levm1, 1)%djnBeg(mm)
627  case (3_inttype)
628  allocate (subface1to1(n1to1)%kdfine(ii), stat=ierr)
629  dfine => subface1to1(n1to1)%kdfine
630  subface1to1(n1to1)%ndk = ii
631  donoroffset = flowdoms(nn, levm1, 1)%dknBeg(mm)
632  end select
633 
634  if (ierr /= 0) &
635  call terminate("createCoarseBlocks", &
636  "Memory allocation failure for jdfine")
637 
638  jj = 1
639  do j = jbeg, jend
640  ii = fact * (jjmap(j) - flowdoms(nn, levm1, 1)%jnBeg(mm))
641  dfine(jj) = ii + donoroffset
642  jj = jj + 1
643  end do
644 
645  ! Determine the fine grid donor indices of the subface.
646  ! k-direction.
647 
648  fact = 1
649  if (l3 < 0) fact = -1
650  ii = kend - kbeg + 1
651 
652  l3 = abs(l3)
653  select case (l3)
654  case (1_inttype)
655  allocate (subface1to1(n1to1)%idfine(ii), stat=ierr)
656  dfine => subface1to1(n1to1)%idfine
657  subface1to1(n1to1)%ndi = ii
658  donoroffset = flowdoms(nn, levm1, 1)%dinBeg(mm)
659  case (2_inttype)
660  allocate (subface1to1(n1to1)%jdfine(ii), stat=ierr)
661  dfine => subface1to1(n1to1)%jdfine
662  subface1to1(n1to1)%ndj = ii
663  donoroffset = flowdoms(nn, levm1, 1)%djnBeg(mm)
664  case (3_inttype)
665  allocate (subface1to1(n1to1)%kdfine(ii), stat=ierr)
666  dfine => subface1to1(n1to1)%kdfine
667  subface1to1(n1to1)%ndk = ii
668  donoroffset = flowdoms(nn, levm1, 1)%dknBeg(mm)
669  end select
670 
671  if (ierr /= 0) &
672  call terminate("createCoarseBlocks", &
673  "Memory allocation failure for kdfine")
674 
675  kk = 1
676  do k = kbeg, kend
677  ii = fact * (kkmap(k) - flowdoms(nn, levm1, 1)%knBeg(mm))
678  dfine(kk) = ii + donoroffset
679  kk = kk + 1
680  end do
681 
682  end if subface_1to1
683 
684  end do subfaces
685 
686  ! Release the local memory allocated inside this loop.
687 
688  deallocate (imap, jmap, kmap, iimap, jjmap, kkmap, stat=ierr)
689  if (ierr /= 0) &
690  call terminate("createCoarseBlocks", &
691  "Deallocation error for imap, iimap, etc.")
692 
693  ! Allocate the memory for the coordinates of all time spectral
694  ! solutions.
695 
696  ii = flowdoms(nn, level, 1)%ie
697  jj = flowdoms(nn, level, 1)%je
698  kk = flowdoms(nn, level, 1)%ke
699 
700  do mm = 1, ntimeintervalsspectral
701  allocate (flowdoms(nn, level, mm)%x(0:ii, 0:jj, 0:kk, 3), stat=ierr)
702  if (ierr /= 0) &
703  call terminate("createCoarseBlocks", &
704  "Memory allocation failure for x")
705  end do
706 
707  ! Copy the scalars such that they are known for all time
708  ! spectral solutions.
709 
710  do mm = 2, ntimeintervalsspectral
711 
712  flowdoms(nn, level, mm)%nx = flowdoms(nn, level, 1)%nx
713  flowdoms(nn, level, mm)%ny = flowdoms(nn, level, 1)%ny
714  flowdoms(nn, level, mm)%nz = flowdoms(nn, level, 1)%nz
715 
716  flowdoms(nn, level, mm)%il = flowdoms(nn, level, 1)%il
717  flowdoms(nn, level, mm)%jl = flowdoms(nn, level, 1)%jl
718  flowdoms(nn, level, mm)%kl = flowdoms(nn, level, 1)%kl
719 
720  flowdoms(nn, level, mm)%ie = flowdoms(nn, level, 1)%ie
721  flowdoms(nn, level, mm)%je = flowdoms(nn, level, 1)%je
722  flowdoms(nn, level, mm)%ke = flowdoms(nn, level, 1)%ke
723 
724  flowdoms(nn, level, mm)%ib = flowdoms(nn, level, 1)%ib
725  flowdoms(nn, level, mm)%jb = flowdoms(nn, level, 1)%jb
726  flowdoms(nn, level, mm)%kb = flowdoms(nn, level, 1)%kb
727 
728  flowdoms(nn, level, mm)%nSubface = flowdoms(nn, level, 1)%nSubface
729  flowdoms(nn, level, mm)%n1to1 = flowdoms(nn, level, 1)%n1to1
730  flowdoms(nn, level, mm)%nBocos = flowdoms(nn, level, 1)%nBocos
731  flowdoms(nn, level, mm)%nViscBocos = flowdoms(nn, level, 1)%nViscBocos
732 
733  flowdoms(nn, level, mm)%iCoarsened = flowdoms(nn, level, 1)%iCoarsened
734  flowdoms(nn, level, mm)%jCoarsened = flowdoms(nn, level, 1)%jCoarsened
735  flowdoms(nn, level, mm)%kCoarsened = flowdoms(nn, level, 1)%kCoarsened
736 
737  flowdoms(nn, level, mm)%blockIsMoving = &
738  flowdoms(nn, level, 1)%blockIsMoving
739  flowdoms(nn, level, mm)%cgnsBlockID = &
740  flowdoms(nn, level, 1)%cgnsBlockID
741 
742  flowdoms(nn, level, mm)%addGridVelocities = &
743  flowdoms(nn, level, 1)%addGridVelocities
744 
745  end do
746 
747  end do domains
748 
749  ! Determine the owned coordinates.
750 
751  call coarseownedcoordinates(level)
752 
753  ! Determine the donor info for the internal block boundaries.
754 
755  call coarsedonorinfo(level)
756 
757  ! Remove possible coarse grid non-matching block boundaries from
758  ! the 1 to 1 list.
759 
760  call checkcoarse1to1(level)
761 
762  ! Release the memory of the coarsening info.
763 
764  do nn = 1, ndom
765  deallocate (coarseinfo(nn)%coarseIs1to1, stat=ierr)
766  if (ierr /= 0) &
767  call terminate("createCoarseBlocks", &
768  "Deallocation error for coarseIs1to1")
769  end do
770 
771  deallocate (coarseinfo, stat=ierr)
772  if (ierr /= 0) &
773  call terminate("createCoarseBlocks", &
774  "Deallocation error for coarseInfo")
775 
776  end subroutine createcoarseblocks
777 
778  ! ==================================================================
779 
780  subroutine coarseownedcoordinates(level)
781  !
782  ! coarseOwnedCoordinates determines from the coarsening info
783  ! the owned coordinates of the coarse grid. This is done in a
784  ! separate routine, because in unsteady moving mesh mode or for
785  ! deforming meshes only new coordinates need to be computed,
786  ! while the connectivity remains the same.
787  !
788  use block
790  implicit none
791  !
792  ! Subroutine arguments.
793  !
794  integer(kind=intType), intent(in) :: level
795  !
796  ! Local variables.
797  !
798  integer(kind=intType) :: i, j, k, ii, jj, kk
799  integer(kind=intType) :: nn, mm, il, jl, kl, levm1
800 
801  logical, dimension(:), pointer :: iCo, jCo, kCo
802 
803  ! Store the finer grid level in levm1
804 
805  levm1 = level - 1
806 
807  ! Loop over the number of blocks on this processor
808 
809  domains: do nn = 1, ndom
810 
811  ! Easier storage of some variables of this block.
812 
813  il = flowdoms(nn, levm1, 1)%il
814  jl = flowdoms(nn, levm1, 1)%jl
815  kl = flowdoms(nn, levm1, 1)%kl
816 
817  ico => flowdoms(nn, levm1, 1)%iCo
818  jco => flowdoms(nn, levm1, 1)%jCo
819  kco => flowdoms(nn, levm1, 1)%kCo
820 
821  ! Loop over the fine grid lines in the three directions and
822  ! determine which should be kept on the coarse grid.
823 
824  kk = 1
825  do k = 1, kl
826  if (kco(k)) then
827  jj = 1
828  do j = 1, jl
829  if (jco(j)) then
830  ii = 1
831  do i = 1, il
832  if (ico(i)) then
833 
834  ! Loop over the spectral solutions and copy the
835  ! coordinates from the fine grid.
836 
837  do mm = 1, ntimeintervalsspectral
838  flowdoms(nn, level, mm)%x(ii, jj, kk, 1) = &
839  flowdoms(nn, levm1, mm)%x(i, j, k, 1)
840  flowdoms(nn, level, mm)%x(ii, jj, kk, 2) = &
841  flowdoms(nn, levm1, mm)%x(i, j, k, 2)
842  flowdoms(nn, level, mm)%x(ii, jj, kk, 3) = &
843  flowdoms(nn, levm1, mm)%x(i, j, k, 3)
844  end do
845 
846  ii = ii + 1
847  end if
848  end do
849  jj = jj + 1
850  end if
851  end do
852  kk = kk + 1
853  end if
854  end do
855 
856  end do domains
857 
858  end subroutine coarseownedcoordinates
859 
860  subroutine update1to1coarse(level, subface)
861  !
862  ! update1to1Coarse determines whether or not the given 1 to 1
863  ! block boundary subface on the fine level is still a 1 to 1
864  ! subface on the coarser grid level.
865  !
866  use block
868  use coarseninginfo
869  use utils, only: terminate
870  implicit none
871  !
872  ! Subroutine arguments.
873  !
874  integer(kind=intType), intent(in) :: level
875  type(coarse1to1subfacetype), intent(in) :: subface
876  !
877  ! Local variables.
878  !
879  integer(kind=intType) :: levm1, nn, mm, ii, jj, kk, ll
880  integer(kind=intType) :: i, j, k
881  integer(kind=intType) :: iBeg, jBeg, kBeg, iEnd, jEnd, kEnd
882 
883  logical :: subfaceFound, idir, jdir, kdir
884 
885  ! Store the finer grid level in levm1
886 
887  levm1 = level - 1
888 
889  ! Store the local block id a bit easier and store the minimum and
890  ! maximum index in the three coordinate directions for the fine
891  ! grid.
892 
893  nn = subface%neighBlock
894 
895  mm = subface%ndi
896  ibeg = min(subface%idfine(1), subface%idfine(mm))
897  iend = max(subface%idfine(1), subface%idfine(mm))
898 
899  mm = subface%ndj
900  jbeg = min(subface%jdfine(1), subface%jdfine(mm))
901  jend = max(subface%jdfine(1), subface%jdfine(mm))
902 
903  mm = subface%ndk
904  kbeg = min(subface%kdfine(1), subface%kdfine(mm))
905  kend = max(subface%kdfine(1), subface%kdfine(mm))
906 
907  ! Find the subface id of the block, which corresponds to the
908  ! given subface of the subroutine header.
909 
910  subfacefound = .false.
911  findsubface: do ii = 1, flowdoms(nn, levm1, 1)%n1to1
912 
913  mm = ii + flowdoms(nn, levm1, 1)%nBocos
914 
915  ! Check if this subface is the correct one.
916 
917  kk = min(flowdoms(nn, levm1, 1)%inBeg(mm), &
918  flowdoms(nn, levm1, 1)%inEnd(mm))
919  ll = max(flowdoms(nn, levm1, 1)%inBeg(mm), &
920  flowdoms(nn, levm1, 1)%inEnd(mm))
921 
922  if (kk == ibeg .and. ll == iend) then
923  kk = min(flowdoms(nn, levm1, 1)%jnBeg(mm), &
924  flowdoms(nn, levm1, 1)%jnEnd(mm))
925  ll = max(flowdoms(nn, levm1, 1)%jnBeg(mm), &
926  flowdoms(nn, levm1, 1)%jnEnd(mm))
927 
928  if (kk == jbeg .and. ll == jend) then
929  kk = min(flowdoms(nn, levm1, 1)%knBeg(mm), &
930  flowdoms(nn, levm1, 1)%knEnd(mm))
931  ll = max(flowdoms(nn, levm1, 1)%knBeg(mm), &
932  flowdoms(nn, levm1, 1)%knEnd(mm))
933 
934  if (kk == kbeg .and. ll == kend) subfacefound = .true.
935  end if
936  end if
937 
938  ! Exit the loop if this is indeed the subface.
939 
940  if (subfacefound) exit
941  end do findsubface
942 
943  ! The subface must be found on the fine grid. Check this.
944 
945  if (.not. subfacefound) &
946  call terminate("update1to1Coarse", &
947  "Invalid fine grid 1 to 1 subface connectivity")
948 
949  ! Check the i-direction of the coarse grid subface to see if it
950  ! is still a 1 to 1 subface. First check the number of grid lines.
951  ! If these are identical, check each coarse grid line.
952 
953  idir = .true.
954  ii = abs(flowdoms(nn, level, 1)%inEnd(mm) &
955  - flowdoms(nn, level, 1)%inBeg(mm)) + 1
956 
957  if (ii == subface%ndi) then
958  do i = 1, subface%ndi
959  ii = subface%idfine(i)
960  if (.not. flowdoms(nn, levm1, 1)%ico(ii)) idir = .false.
961  end do
962  else
963  idir = .false.
964  end if
965 
966  ! Check the j-direction of the coarse grid subface to see if it
967  ! is still a 1 to 1 subface. First check the number of grid lines.
968  ! If these are identical, check each coarse grid line.
969 
970  jdir = .true.
971  jj = abs(flowdoms(nn, level, 1)%jnEnd(mm) &
972  - flowdoms(nn, level, 1)%jnBeg(mm)) + 1
973 
974  if (jj == subface%ndj) then
975  do j = 1, subface%ndj
976  jj = subface%jdfine(j)
977  if (.not. flowdoms(nn, levm1, 1)%jco(jj)) jdir = .false.
978  end do
979  else
980  jdir = .false.
981  end if
982 
983  ! Check the k-direction of the coarse grid subface to see if it
984  ! is still a 1 to 1 subface. First check the number of grid lines.
985  ! If these are identical, check each coarse grid line.
986 
987  kdir = .true.
988  kk = abs(flowdoms(nn, level, 1)%knEnd(mm) &
989  - flowdoms(nn, level, 1)%knBeg(mm)) + 1
990 
991  if (kk == subface%ndk) then
992  do k = 1, subface%ndk
993  kk = subface%kdfine(k)
994  if (.not. flowdoms(nn, levm1, 1)%kco(kk)) kdir = .false.
995  end do
996  else
997  kdir = .false.
998  end if
999 
1000  ! Set coarseIs1to1 to .true. if the subface on the coarse grid
1001  ! is still a 1 to 1 subface. Otherwise set it to .false.
1002 
1003  ii = mm - flowdoms(nn, levm1, 1)%nBocos
1004  if (idir .and. jdir .and. kdir) then
1005  coarseinfo(nn)%coarseIs1to1(ii) = .true.
1006  else
1007  coarseinfo(nn)%coarseIs1to1(ii) = .false.
1008  end if
1009 
1010  ! Store the donor range. It is possible that the lower and
1011  ! upper boundary must be reversed.
1012 
1013  if (flowdoms(nn, levm1, 1)%dinBeg(mm) < &
1014  flowdoms(nn, levm1, 1)%dinEnd(mm)) then
1015  flowdoms(nn, level, 1)%dinBeg(mm) = subface%iBeg
1016  flowdoms(nn, level, 1)%dinEnd(mm) = subface%iEnd
1017  else
1018  flowdoms(nn, level, 1)%dinBeg(mm) = subface%iEnd
1019  flowdoms(nn, level, 1)%dinEnd(mm) = subface%iBeg
1020  end if
1021 
1022  if (flowdoms(nn, levm1, 1)%djnBeg(mm) < &
1023  flowdoms(nn, levm1, 1)%djnEnd(mm)) then
1024  flowdoms(nn, level, 1)%djnBeg(mm) = subface%jBeg
1025  flowdoms(nn, level, 1)%djnEnd(mm) = subface%jEnd
1026  else
1027  flowdoms(nn, level, 1)%djnBeg(mm) = subface%jEnd
1028  flowdoms(nn, level, 1)%djnEnd(mm) = subface%jBeg
1029  end if
1030 
1031  if (flowdoms(nn, levm1, 1)%dknBeg(mm) < &
1032  flowdoms(nn, levm1, 1)%dknEnd(mm)) then
1033  flowdoms(nn, level, 1)%dknBeg(mm) = subface%kBeg
1034  flowdoms(nn, level, 1)%dknEnd(mm) = subface%kEnd
1035  else
1036  flowdoms(nn, level, 1)%dknBeg(mm) = subface%kEnd
1037  flowdoms(nn, level, 1)%dknEnd(mm) = subface%kBeg
1038  end if
1039 
1040  end subroutine update1to1coarse
1041  subroutine coarsedonorinfo(level)
1042  !
1043  ! coarseDonorInfo creates the donor info for the internal
1044  ! block boundaries on the given coarse grid level.
1045  !
1046  use communication
1047  use coarse1to1subface
1048  use utils, only: terminate
1049  implicit none
1050  !
1051  ! Subroutine arguments.
1052  !
1053  integer(kind=intType), intent(in) :: level
1054  !
1055  ! Local variables.
1056  !
1057  integer :: proc, size, ierr
1058 
1059  integer, dimension(mpi_status_size) :: mpiStatus
1060  integer, dimension(nProc) :: sizeMessage
1061 
1062  integer(kind=intType) :: nn, mm, ii, i
1063  integer(kind=intType) :: nMessageReceive, nMessageSend
1064 
1065  integer(kind=intType), dimension(0:nProc) :: size2proc
1066  integer(kind=intType), dimension(nProc) :: tmp
1067 
1068  integer(kind=intType), dimension(:), allocatable :: sendBuf
1069  integer(kind=intType), dimension(:), allocatable, target :: recvBuf
1070 
1071  type(coarse1to1subfacetype) :: tmpSubface
1072 
1073  ! Determine the number of integers sent to every processor.
1074 
1075  size2proc = 0
1076  do nn = 1, nsubface1to1
1077  mm = subface1to1(nn)%neighProc + 1 ! Proc ID's start at 0.
1078 
1079  size2proc(mm) = size2proc(mm) + 13 &
1080  + subface1to1(nn)%iEnd - subface1to1(nn)%iBeg &
1081  + subface1to1(nn)%jEnd - subface1to1(nn)%jBeg &
1082  + subface1to1(nn)%kEnd - subface1to1(nn)%kBeg
1083  end do
1084 
1085  ! No message is sent to myself, so set the corresponding entry in
1086  ! size2proc to 0. The processor id's start at 0, which explains
1087  ! the +1.
1088 
1089  size2proc(myid + 1) = 0
1090 
1091  ! Determine the number of messages i have to sent. Store in tmp
1092  ! a 0 or a 1, depending whether or not a message must be sent
1093  ! to the corresponding processor. This info is needed to determine
1094  ! the number of messages i will receive.
1095 
1096  nmessagesend = 0
1097  do nn = 1, nproc
1098  if (size2proc(nn) > 0) then
1099  nmessagesend = nmessagesend + 1
1100  tmp(nn) = 1
1101  else
1102  tmp(nn) = 0
1103  end if
1104  end do
1105 
1106  ! Determine the number of messages i will receive.
1107 
1108  sizemessage = 1
1109  call mpi_reduce_scatter(tmp, nmessagereceive, sizemessage, &
1110  adflow_integer, mpi_sum, adflow_comm_world, &
1111  ierr)
1112 
1113  ! Put size2proc in cumulative storage format and store the
1114  ! starting entry in tmp, which is used as a counter.
1115 
1116  do nn = 1, nproc
1117  size2proc(nn) = size2proc(nn) + size2proc(nn - 1)
1118  tmp(nn) = size2proc(nn - 1)
1119  end do
1120 
1121  ! Allocate the memory for the send buffer.
1122 
1123  allocate (sendbuf(size2proc(nproc)), stat=ierr)
1124  if (ierr /= 0) &
1125  call terminate("coarseDonorInfo", &
1126  "Memory allocation failure for sendBuf")
1127 
1128  ! Loop over the number of 1 to 1 subfaces to fill the send buffer.
1129 
1130  unownedsubfaces: do nn = 1, nsubface1to1
1131 
1132  ! Only info that must be sent to other processors must be stored.
1133 
1134  if (subface1to1(nn)%neighProc /= myid) then
1135 
1136  ! Store the entry in tmp in mm.
1137  ! Note that the proc id's start at 0.
1138 
1139  mm = subface1to1(nn)%neighProc + 1
1140 
1141  ! Store the block id of the donor and the coarse grid range of
1142  ! the current subface in the send buffer. Update the counter
1143  ! tmp(mm) accordingly.
1144 
1145  tmp(mm) = tmp(mm) + 1
1146  sendbuf(tmp(mm)) = subface1to1(nn)%neighBlock
1147 
1148  tmp(mm) = tmp(mm) + 1
1149  sendbuf(tmp(mm)) = subface1to1(nn)%iBeg
1150 
1151  tmp(mm) = tmp(mm) + 1
1152  sendbuf(tmp(mm)) = subface1to1(nn)%jBeg
1153 
1154  tmp(mm) = tmp(mm) + 1
1155  sendbuf(tmp(mm)) = subface1to1(nn)%kBeg
1156 
1157  tmp(mm) = tmp(mm) + 1
1158  sendbuf(tmp(mm)) = subface1to1(nn)%iEnd
1159 
1160  tmp(mm) = tmp(mm) + 1
1161  sendbuf(tmp(mm)) = subface1to1(nn)%jEnd
1162 
1163  tmp(mm) = tmp(mm) + 1
1164  sendbuf(tmp(mm)) = subface1to1(nn)%kEnd
1165 
1166  ! Store the number of points in the three coordinate
1167  ! directions of the coarse grid donor face.
1168 
1169  tmp(mm) = tmp(mm) + 1
1170  sendbuf(tmp(mm)) = subface1to1(nn)%ndi
1171 
1172  tmp(mm) = tmp(mm) + 1
1173  sendbuf(tmp(mm)) = subface1to1(nn)%ndj
1174 
1175  tmp(mm) = tmp(mm) + 1
1176  sendbuf(tmp(mm)) = subface1to1(nn)%ndk
1177 
1178  ! Store the i-donor indices of the fine grid in sendBuf.
1179 
1180  do i = 1, subface1to1(nn)%ndi
1181  tmp(mm) = tmp(mm) + 1
1182  sendbuf(tmp(mm)) = subface1to1(nn)%idfine(i)
1183  end do
1184 
1185  ! Store the j-donor indices of the fine grid in sendBuf.
1186 
1187  do i = 1, subface1to1(nn)%ndj
1188  tmp(mm) = tmp(mm) + 1
1189  sendbuf(tmp(mm)) = subface1to1(nn)%jdfine(i)
1190  end do
1191 
1192  ! Store the k-donor indices of the fine grid in sendBuf.
1193 
1194  do i = 1, subface1to1(nn)%ndk
1195  tmp(mm) = tmp(mm) + 1
1196  sendbuf(tmp(mm)) = subface1to1(nn)%kdfine(i)
1197  end do
1198 
1199  end if
1200 
1201  end do unownedsubfaces
1202 
1203  ! Send the data i have to send.
1204 
1205  mm = 0
1206  sends: do nn = 1, nproc
1207 
1208  if (size2proc(nn) > size2proc(nn - 1)) then
1209 
1210  ! Update mm and store the processor id, the size of the
1211  ! message and the starting index in sendbuf.
1212 
1213  mm = mm + 1
1214  proc = nn - 1
1215  size = size2proc(nn) - size2proc(nn - 1)
1216  ii = size2proc(nn - 1) + 1
1217 
1218  ! Send the message.
1219 
1220  call mpi_isend(sendbuf(ii), size, adflow_integer, proc, proc, &
1221  adflow_comm_world, sendrequests(mm), ierr)
1222 
1223  end if
1224 
1225  end do sends
1226 
1227  ! Loop over the number of 1 to 1 subfaces stored on this processor
1228  ! and update the locally stored info.
1229 
1230  do nn = 1, nsubface1to1
1231  if (subface1to1(nn)%neighProc == myid) &
1232  call update1to1coarse(level, subface1to1(nn))
1233  end do
1234 
1235  ! Release the memory of subface1to1.
1236 
1237  do nn = 1, nsubface1to1
1238  deallocate (subface1to1(nn)%idfine, subface1to1(nn)%jdfine, &
1239  subface1to1(nn)%kdfine, stat=ierr)
1240  if (ierr /= 0) &
1241  call terminate("coarseDonorInfo", &
1242  "Deallocation error for idfine, etc")
1243  end do
1244 
1245  deallocate (subface1to1, stat=ierr)
1246  if (ierr /= 0) &
1247  call terminate("coarseDonorInfo", &
1248  "Deallocation error for subface1to1")
1249 
1250  ! Loop over the number of messages i must receive to determine
1251  ! info from externally stored neighboring subfaces.
1252 
1253  receives: do nn = 1, nmessagereceive
1254 
1255  ! Block until a message arrives.
1256 
1257  call mpi_probe(mpi_any_source, myid, adflow_comm_world, &
1258  mpistatus, ierr)
1259 
1260  ! Find the source and size of the message.
1261 
1262  proc = mpistatus(mpi_source)
1263  call mpi_get_count(mpistatus, adflow_integer, size, ierr)
1264 
1265  ! Check in debug mode that the incoming message is of
1266  ! correct size.
1267 
1268  if (debug) then
1269  if (size == mpi_undefined) &
1270  call terminate("coarseDonorInfo", &
1271  "Unexpected size of message")
1272  end if
1273 
1274  ! Allocate the memory for the receive buffer.
1275 
1276  allocate (recvbuf(size), stat=ierr)
1277  if (ierr /= 0) &
1278  call terminate("coarseDonorInfo", &
1279  "Memory allocation failure for recvBuf")
1280 
1281  ! Receive the messsage. As it has already arrived a blocking
1282  ! receive can be used.
1283 
1284  call mpi_recv(recvbuf, size, adflow_integer, proc, myid, &
1285  adflow_comm_world, mpistatus, ierr)
1286 
1287  ! Loop to extract the 1 to 1 subface info from the buffer.
1288 
1289  ii = 1
1290  extractsubface: do
1291  ! Exit the loop if no more info is present.
1292 
1293  if (ii > size) exit
1294 
1295  ! Store the info of this subface in tmpSubface.
1296 
1297  tmpsubface%neighProc = proc
1298  tmpsubface%neighBlock = recvbuf(ii); ii = ii + 1
1299 
1300  tmpsubface%iBeg = recvbuf(ii); ii = ii + 1
1301  tmpsubface%jBeg = recvbuf(ii); ii = ii + 1
1302  tmpsubface%kBeg = recvbuf(ii); ii = ii + 1
1303 
1304  tmpsubface%iEnd = recvbuf(ii); ii = ii + 1
1305  tmpsubface%jEnd = recvbuf(ii); ii = ii + 1
1306  tmpsubface%kEnd = recvbuf(ii); ii = ii + 1
1307 
1308  tmpsubface%ndi = recvbuf(ii); ii = ii + 1
1309  tmpsubface%ndj = recvbuf(ii); ii = ii + 1
1310  tmpsubface%ndk = recvbuf(ii); ii = ii + 1
1311 
1312  mm = ii + tmpsubface%ndi
1313  tmpsubface%idfine => recvbuf(ii:(mm - 1))
1314  ii = mm
1315 
1316  mm = ii + tmpsubface%ndj
1317  tmpsubface%jdfine => recvbuf(ii:(mm - 1))
1318  ii = mm
1319 
1320  mm = ii + tmpsubface%ndk
1321  tmpsubface%kdfine => recvbuf(ii:(mm - 1))
1322  ii = mm
1323 
1324  ! Update the 1 to 1 subface info.
1325 
1326  call update1to1coarse(level, tmpsubface)
1327 
1328  end do extractsubface
1329 
1330  ! Release the memory of the receive buffer.
1331 
1332  deallocate (recvbuf, stat=ierr)
1333  if (ierr /= 0) &
1334  call terminate("coarseDonorInfo", &
1335  "Deallocation failure for recvBuf")
1336 
1337  end do receives
1338 
1339  ! Complete the nonblocking sends.
1340 
1341  size = nmessagesend
1342  do nn = 1, nmessagesend
1343  call mpi_waitany(size, sendrequests, proc, mpistatus, ierr)
1344  end do
1345 
1346  ! Release the memory of the send buffer.
1347 
1348  deallocate (sendbuf, stat=ierr)
1349  if (ierr /= 0) &
1350  call terminate("coarseDonorInfo", &
1351  "Deallocation failure for sendBuf")
1352 
1353  ! Synchronize the processors, because wild cards have been used.
1354 
1355  call mpi_barrier(adflow_comm_world, ierr)
1356 
1357  end subroutine coarsedonorinfo
1358 
1359  subroutine checkcoarse1to1(level)
1360  !
1361  ! checkCoarse1to1 removes the nonmatching block boundaries
1362  ! from the list of 1 to 1 matching ones. They are in there,
1363  ! because they are 1 to 1 matching on the finer grids.
1364  !
1365  use constants
1366  use cgnsgrid
1367  use block
1368  use inputtimespectral
1369  use coarseninginfo
1370  use utils, only: terminate
1371  use commonformats, only: strings
1372  implicit none
1373  !
1374  ! Subroutine arguments.
1375  !
1376  integer(kind=intType), intent(in) :: level
1377  !
1378  ! Local variables.
1379  !
1380  character(len=maxStringLen) :: errorMessage
1381 
1382  integer(kind=intType) :: i, nn, mm, kk, ll
1383 
1384  ! Loop over the number of domains.
1385 
1386  domains: do nn = 1, ndom
1387 
1388  ! Loop over the number of 1 to 1 subfaces. As this number can
1389  ! change during the loop, the control statement is done via
1390  ! an exit.
1391 
1392  i = 1
1393  n1to1: do
1394 
1395  ! Exit the loop if the counter is larger than the number
1396  ! of 1 to 1 subfaces.
1397 
1398  if (i > flowdoms(nn, level, 1)%n1to1) exit
1399 
1400  ! Add the offset to i to store the correct place in the
1401  ! subface info.
1402 
1403  mm = i + flowdoms(nn, level, 1)%nBocos
1404 
1405  ! Test if this is still a 1 to 1 subface on the coarse grid.
1406 
1407  is1to1: if (coarseinfo(nn)%coarseIs1to1(i)) then
1408 
1409  ! Subface is still a 1 to 1 subface on the coarse grid.
1410  ! Only the counter i must be updated.
1411 
1412  i = i + 1
1413 
1414  else is1to1
1415 
1416  ! Due to the coarsening the subface is not a 1 to 1 block
1417  ! boundary anymore. Swap the current entry with the entry
1418  ! kk.
1419 
1420  kk = flowdoms(nn, level, 1)%nBocos + flowdoms(nn, level, 1)%n1to1
1421 
1422  ! The range info.
1423 
1424  ll = flowdoms(nn, level, 1)%inBeg(mm)
1425  flowdoms(nn, level, 1)%inBeg(mm) = flowdoms(nn, level, 1)%inBeg(kk)
1426  flowdoms(nn, level, 1)%inBeg(kk) = ll
1427 
1428  ll = flowdoms(nn, level, 1)%jnBeg(mm)
1429  flowdoms(nn, level, 1)%jnBeg(mm) = flowdoms(nn, level, 1)%jnBeg(kk)
1430  flowdoms(nn, level, 1)%jnBeg(kk) = ll
1431 
1432  ll = flowdoms(nn, level, 1)%knBeg(mm)
1433  flowdoms(nn, level, 1)%knBeg(mm) = flowdoms(nn, level, 1)%knBeg(kk)
1434  flowdoms(nn, level, 1)%knBeg(kk) = ll
1435 
1436  ll = flowdoms(nn, level, 1)%inEnd(mm)
1437  flowdoms(nn, level, 1)%inEnd(mm) = flowdoms(nn, level, 1)%inEnd(kk)
1438  flowdoms(nn, level, 1)%inEnd(kk) = ll
1439 
1440  ll = flowdoms(nn, level, 1)%jnEnd(mm)
1441  flowdoms(nn, level, 1)%jnEnd(mm) = flowdoms(nn, level, 1)%jnEnd(kk)
1442  flowdoms(nn, level, 1)%jnEnd(kk) = ll
1443 
1444  ll = flowdoms(nn, level, 1)%knEnd(mm)
1445  flowdoms(nn, level, 1)%knEnd(mm) = flowdoms(nn, level, 1)%knEnd(kk)
1446  flowdoms(nn, level, 1)%knEnd(kk) = ll
1447 
1448  ! The donor info.
1449 
1450  ll = flowdoms(nn, level, 1)%dinBeg(mm)
1451  flowdoms(nn, level, 1)%dinBeg(mm) = flowdoms(nn, level, 1)%dinBeg(kk)
1452  flowdoms(nn, level, 1)%dinBeg(kk) = ll
1453 
1454  ll = flowdoms(nn, level, 1)%djnBeg(mm)
1455  flowdoms(nn, level, 1)%djnBeg(mm) = flowdoms(nn, level, 1)%djnBeg(kk)
1456  flowdoms(nn, level, 1)%djnBeg(kk) = ll
1457 
1458  ll = flowdoms(nn, level, 1)%dknBeg(mm)
1459  flowdoms(nn, level, 1)%dknBeg(mm) = flowdoms(nn, level, 1)%dknBeg(kk)
1460  flowdoms(nn, level, 1)%dknBeg(kk) = ll
1461 
1462  ll = flowdoms(nn, level, 1)%dinEnd(mm)
1463  flowdoms(nn, level, 1)%dinEnd(mm) = flowdoms(nn, level, 1)%dinEnd(kk)
1464  flowdoms(nn, level, 1)%dinEnd(kk) = ll
1465 
1466  ll = flowdoms(nn, level, 1)%djnEnd(mm)
1467  flowdoms(nn, level, 1)%djnEnd(mm) = flowdoms(nn, level, 1)%djnEnd(kk)
1468  flowdoms(nn, level, 1)%djnEnd(kk) = ll
1469 
1470  ll = flowdoms(nn, level, 1)%dknEnd(mm)
1471  flowdoms(nn, level, 1)%dknEnd(mm) = flowdoms(nn, level, 1)%dknEnd(kk)
1472  flowdoms(nn, level, 1)%dknEnd(kk) = ll
1473 
1474  ! The transformation matrix.
1475 
1476  ll = flowdoms(nn, level, 1)%l1(mm)
1477  flowdoms(nn, level, 1)%l1(mm) = flowdoms(nn, level, 1)%l1(kk)
1478  flowdoms(nn, level, 1)%l1(kk) = ll
1479 
1480  ll = flowdoms(nn, level, 1)%l2(mm)
1481  flowdoms(nn, level, 1)%l2(mm) = flowdoms(nn, level, 1)%l2(kk)
1482  flowdoms(nn, level, 1)%l2(kk) = ll
1483 
1484  ll = flowdoms(nn, level, 1)%l3(mm)
1485  flowdoms(nn, level, 1)%l3(mm) = flowdoms(nn, level, 1)%l3(kk)
1486  flowdoms(nn, level, 1)%l3(kk) = ll
1487 
1488  ! The rest of the subface info.
1489 
1490  ll = flowdoms(nn, level, 1)%BCType(mm)
1491  flowdoms(nn, level, 1)%BCType(mm) = &
1492  flowdoms(nn, level, 1)%BCType(kk)
1493  flowdoms(nn, level, 1)%BCType(kk) = ll
1494 
1495  ll = flowdoms(nn, level, 1)%BCFaceID(mm)
1496  flowdoms(nn, level, 1)%BCFaceID(mm) = &
1497  flowdoms(nn, level, 1)%BCFaceID(kk)
1498  flowdoms(nn, level, 1)%BCFaceID(kk) = ll
1499 
1500  ll = flowdoms(nn, level, 1)%neighBlock(mm)
1501  flowdoms(nn, level, 1)%neighBlock(mm) = &
1502  flowdoms(nn, level, 1)%neighBlock(kk)
1503  flowdoms(nn, level, 1)%neighBlock(kk) = ll
1504 
1505  ll = flowdoms(nn, level, 1)%neighProc(mm)
1506  flowdoms(nn, level, 1)%neighProc(mm) = &
1507  flowdoms(nn, level, 1)%neighProc(kk)
1508  flowdoms(nn, level, 1)%neighProc(kk) = ll
1509 
1510  ll = flowdoms(nn, level, 1)%groupNum(mm)
1511  flowdoms(nn, level, 1)%groupNum(mm) = &
1512  flowdoms(nn, level, 1)%groupNum(kk)
1513  flowdoms(nn, level, 1)%groupNum(kk) = ll
1514 
1515  ll = flowdoms(nn, level, 1)%cgnsSubface(mm)
1516  flowdoms(nn, level, 1)%cgnsSubface(mm) = &
1517  flowdoms(nn, level, 1)%cgnsSubface(kk)
1518  flowdoms(nn, level, 1)%cgnsSubface(kk) = ll
1519 
1520  ! Decrease the number of 1 to 1 block boundaries. Note the
1521  ! counter i should not be updated.
1522 
1523  flowdoms(nn, level, 1)%n1to1 = flowdoms(nn, level, 1)%n1to1 - 1
1524 
1525  ll = flowdoms(nn, 1, 1)%cgnsBlockID
1526  write (errormessage, strings) "Non-matching block-to-block face on zone ", cgnsdoms(ll)%zoneName, &
1527  "...Support not implemented yet."
1528  call terminate("checkCoarse1to1", errormessage)
1529 
1530  end if is1to1
1531 
1532  end do n1to1
1533 
1534  ! Copy the number of internal subfaces for the rest of
1535  ! the spectral solutions.
1536 
1537  do i = 2, ntimeintervalsspectral
1538  flowdoms(nn, level, i)%n1to1 = flowdoms(nn, level, 1)%n1to1
1539  end do
1540 
1541  end do domains
1542 
1543  end subroutine checkcoarse1to1
1544 
1545 end module coarseutils
Definition: block.F90:1
integer(kind=inttype) ndom
Definition: block.F90:761
type(blocktype), dimension(:, :, :), allocatable, target flowdoms
Definition: block.F90:771
type(cgnsblockinfotype), dimension(:), allocatable cgnsdoms
Definition: cgnsGrid.F90:495
type(coarse1to1subfacetype), dimension(:), allocatable subface1to1
integer(kind=inttype) nsubface1to1
type(coarseninginfotype), dimension(:), allocatable coarseinfo
subroutine createcoarseblocks(level)
Definition: coarseUtils.F90:6
subroutine update1to1coarse(level, subface)
subroutine checkcoarse1to1(level)
subroutine coarsedonorinfo(level)
subroutine coarseownedcoordinates(level)
character(len=maxstringlen) strings
integer, dimension(:), allocatable sendrequests
integer adflow_comm_world
integer(kind=portype), parameter regular
Definition: constants.F90:229
real(kind=realtype), parameter one
Definition: constants.F90:72
real(kind=realtype), parameter half
Definition: constants.F90:80
integer(kind=portype), parameter rightstarted
Definition: constants.F90:229
integer(kind=portype), parameter leftstarted
Definition: constants.F90:229
integer(kind=inttype) ntimeintervalsspectral
Definition: inputParam.F90:645
Definition: utils.F90:1
subroutine terminate(routineName, errorMessage)
Definition: utils.F90:502