Line data Source code
1 : !> \file mo_append.f90
2 : !> \brief \copybrief mo_append
3 : !> \details \copydetails mo_append
4 :
5 : !> \brief Append values on existing arrays.
6 : !> \details Provides routines to append (rows) and paste (columns) scalars, vectors,
7 : !! and matrixes onto existing arrays.
8 : !> \changelog
9 : !! - Juliane Mai, Aug 2012
10 : !! - Juliane Mai, Aug 2012
11 : !! - character append & paste
12 : !! - Matthias Cuntz, Jan 2013
13 : !! - removed 256 character restriction
14 : !! - Matthias Cuntz, Feb 2013
15 : !! - logical append and paste
16 : !> \author Juliane Mai
17 : !> \date Aug 2012
18 : !> \copyright Copyright 2005-\today, the CHS Developers, Sabine Attinger: All rights reserved.
19 : !! FORCES is released under the LGPLv3+ license \license_note
20 : MODULE mo_append
21 :
22 : USE mo_kind, only: i4, i8, sp, dp
23 :
24 : IMPLICIT NONE
25 :
26 : PUBLIC :: append !< Returns input1 appended (on first dim) with input2. (like Unix cat)
27 : PUBLIC :: paste !< Returns input1 pasted (on last dim) with input2. (like Unix paste)
28 : PUBLIC :: add_nodata_slice !< Returns input1 extended (on last dim) with nodata slice.
29 :
30 : ! ------------------------------------------------------------------
31 :
32 : !> \brief Append (rows) scalars, vectors, and matrixes onto existing array.
33 :
34 : !> \details Appends one input to the rows of another, i.e. append
35 : !! on the first dimension. Append variants for 3d and 4d also exists.\n
36 : !! The input might be a scalar, a vector or a matrix.\n
37 : !! Possibilities are:\n
38 : !! 1. append scalar to vector\n
39 : !! 2. append vector to vector\n
40 : !! 3. append matrix to matrix
41 : !!
42 : !> \b Example
43 : !!
44 : !! Append input1 with input2
45 : !! \code{.f90}
46 : !! input1 = (/ 1.0_dp , 2.0_dp /)
47 : !! input2 = 3.0_dp
48 : !!
49 : !! call append(input1, input2)
50 : !! --> input1 = (/ 1.0_dp , 2.0_dp, 3.0_dp /)
51 : !! \endcode
52 : !!
53 : !! See also test folder for a detailed example, "test/test_mo_append/".
54 :
55 : !> \param[in] "input2"
56 : !! Values to append. Can be `integer(i4/i8)`, `real(sp/dp)`, `logical`, or `character(len=*)`
57 : !! and also scalar, `dimension(:)`, or `dimension(:,:)`.\n
58 : !! Includes 3d version, where the append is always done along
59 : !! the first dimension.\n
60 : !! 4d version is available for `real(dp)` and `logical`, where
61 : !! append can be done in any component of the dimension.\n
62 : !! If not scalar then the columns have to agree with `input1`.
63 : !> \param[inout] "`allocatable` :: input1"
64 : !! Array to be appended to. Can be `integer(i4/i8)`, `real(sp/dp)`, `logical`,
65 : !! or `character(len=*)`. Must be `dimension(:)` or `dimension(:,:)`, and `allocatable`.\n
66 : !! If `input2` is not scalar then it must be `size(input1,2) = size(input2,2)`.
67 : !> \param[in] "optional :: fill_value"
68 : !! Filler value if number of columns do not coincide. Then the largest column will
69 : !! be taken and empty elements from the generated column will be filled by this value.
70 : !> \param[in] "integer(i4), optional :: idim"
71 : !! Applicable when appending matrix-matrix, 3d, and 4d for `logical` and `real(dp)`.
72 : !! Determines along which dimension `input2` is appended to.
73 :
74 : !> \note
75 : !! Size of `input1` and `input2` have to fit together,
76 : !! i.e. number of columns input1 = number of columns input2.
77 :
78 : !> \author Juliane Mai
79 : !> \date Aug 2012
80 :
81 : !> \author Matthias Cuntz
82 : !> \date Jan 2013
83 : !! - removed 256 character restriction.
84 : !> \date Feb 2013
85 : !! - logical append and paste.
86 :
87 : !> \author Matthias Zink
88 : !> \date Feb 2015
89 : !! - added optional 'fill_value' for logical append
90 :
91 : !> \author Stephan Thober
92 : !> \date Jan 2017
93 : !! - added 3d version for append
94 : !> \date Jul 2018
95 : !! - added optional iDim argument for arrays larger equal 2d
96 :
97 : !> \author Arya Prasetya
98 : !> \date Nov 2021
99 : !! - included i4_3d in append interface
100 :
101 : INTERFACE append
102 : MODULE PROCEDURE append_i4_v_s, append_i4_v_v, append_i4_m_m, &
103 : append_i4_3d, append_i8_v_s, append_i8_v_v, &
104 : append_i8_m_m, append_i8_3d, &
105 : append_sp_v_s, append_sp_v_v, append_sp_m_m, &
106 : append_sp_3d, &
107 : append_dp_v_s, append_dp_v_v, append_dp_m_m, &
108 : append_dp_3d, append_dp_4d, &
109 : append_char_v_s, append_char_v_v, append_char_m_m, &
110 : append_char_3d, &
111 : append_lgt_v_s, append_lgt_v_v, append_lgt_m_m, &
112 : append_lgt_3d, append_lgt_4d
113 : END INTERFACE append
114 :
115 : ! ------------------------------------------------------------------
116 :
117 : !> \brief Paste (columns) scalars, vectors, and matrixes onto existing array.
118 :
119 : !> \details Pastes one input to the columns of another, i.e. append
120 : !! on the second dimension. Paste variants for 3d and 4d also exists.\n
121 : !! The input might be a scalar, a vector or a matrix.\n
122 : !! Possibilities are:\n
123 : !! 1. paste scalar to one-line matrix\n
124 : !! 2. paste vector to a matrix\n
125 : !! 3. paste matrix to matrix
126 : !!
127 : !! \b Example
128 : !!
129 : !! Paste input2 to input1.
130 : !! \code{.f90}
131 : !! input1 = (/ 1.0_dp , 2.0_dp /)
132 : !! input2 = (/ 3.0_dp , 4.0_dp /)
133 : !!
134 : !! call paste(input1, input2)
135 : !! --> input1(1,:) = (/ 1.0_dp , 3.0_dp /)
136 : !! input1(2,:) = (/ 2.0_dp , 4.0_dp /)
137 : !! \endcode
138 : !! See also test folder for a detailed example, "test/test_mo_append".
139 :
140 : !> \param[in] "input2"
141 : !! Values to paste. Can be `integer(i4/i8)`, `real(sp/dp)`, `logical`, or `character(len=*)`
142 : !! and also scalar, `dimension(:)`, or `dimension(:,:)`.\n
143 : !! Includes 3d and 4d version for `real(dp)` and `integer(i4)`, where paste is always done along
144 : !! the last dimension.\n
145 : !! If not scalar then the rows have to agree with `input1`.\n
146 : !> \param[inout] "allocatable :: input1"
147 : !! Array to be pasted to. Can be `integer(i4/i8)`, `real(sp/dp)`,
148 : !! or `character(len=*)`. Must be `dimension(:)`, `dimension(:,:)`, `dimension(:,:,:)`, or
149 : !! `dimension(:,:,:,:)` and `allocatable`.\n
150 : !! If `input2` is not scalar then it must be `size(input1,1) = size(input2,1)`.
151 : !> \param[in] "optional :: fill_value"
152 : !! Filler value if number of rows do not coincide. Then the largest row will
153 : !! be taken and empty elements from the generated rows will be filled by this value.
154 :
155 : !> \note
156 : !! Size of input1 and input2 have to fit together,
157 : !! i.e. number of rows input1 = number of rows input2
158 :
159 : !> \author Juliane Mai
160 : !> \date Aug 2012
161 :
162 : !> \author Matthias Cuntz
163 : !> \date Jan 2013
164 : !! - removed 256 character restriction
165 : !> \date Feb 2013
166 : !! - logical append and paste
167 :
168 : INTERFACE paste
169 : MODULE PROCEDURE paste_i4_m_s, paste_i4_m_v, paste_i4_m_m, &
170 : paste_i8_m_s, paste_i8_m_v, paste_i8_m_m, &
171 : paste_sp_m_s, paste_sp_m_v, paste_sp_m_m, &
172 : paste_dp_m_s, paste_dp_m_v, paste_dp_m_m, &
173 : paste_char_m_s, paste_char_m_v, paste_char_m_m, &
174 : paste_lgt_m_s, paste_lgt_m_v, paste_lgt_m_m, &
175 : paste_dp_3d, paste_dp_4d, paste_i4_3d, paste_i4_4d
176 :
177 : END INTERFACE paste
178 :
179 : ! ------------------------------------------------------------------
180 :
181 : !> \brief Paste a matrix of ones times a value onto an existing matrix.
182 :
183 : !> \details Pastes a matrix of uniform element and predefined number of columns
184 : !! to the columns of another matrix.\n
185 : !!
186 : !! \b Example
187 : !!
188 : !! Add 2 columns of elements 5.0_dp.
189 : !! \code{.f90}
190 : !! input1(1,:) = (/ 1.0_dp , 2.0_dp /)
191 : !! input1(2,:) = (/ 2.0_dp , 4.0_dp /)
192 : !!
193 : !! call add_nodata_slice(input1, 2, 5.0_dp)
194 : !! --> input1(1,:) = (/ 1.0_dp , 2.0_dp /)
195 : !! input1(2,:) = (/ 2.0_dp , 4.0_dp /)
196 : !! input1(3,:) = (/ 5.0_dp , 5.0_dp /)
197 : !! input1(4,:) = (/ 5.0_dp , 5.0_dp /)
198 : !! \endcode
199 : !!
200 : !! See also test folder for a detailed example, "test/test_mo_append/".
201 :
202 : !> \param[in,out] "allocatable :: matrix"
203 : !! Matrix to be pasted to. Can be `integer(i4)` or `real(dp)`.
204 : !! Must be `dimension(:,:)`, `dimension(:,:,:)`, or `dimension(:,:,:,:)` and `allocatable`.\n
205 : !> \param[in] "integer(i4) :: nAdd"
206 : !! Number of column to paste.\n
207 : !> \param[in] "noDataValue"
208 : !! Value of elements of the matrix. Can be `integer(i4)` or `real(dp)`.
209 :
210 : interface add_nodata_slice
211 : module procedure add_nodata_slice_dp_2d, add_nodata_slice_dp_3d, add_nodata_slice_dp_4d, &
212 : add_nodata_slice_i4_2d, add_nodata_slice_i4_3d, add_nodata_slice_i4_4d
213 : end interface add_nodata_slice
214 :
215 : ! ------------------------------------------------------------------
216 :
217 : PRIVATE
218 :
219 : ! ------------------------------------------------------------------
220 :
221 : CONTAINS
222 :
223 : ! ------------------------------------------------------------------
224 :
225 3 : SUBROUTINE append_i4_v_s(vec1, sca2)
226 :
227 : implicit none
228 :
229 : integer(i4), dimension(:), allocatable, intent(inout) :: vec1
230 : integer(i4), intent(in) :: sca2
231 :
232 : ! local variables
233 : integer(i4) :: n1, n2
234 3 : integer(i4), dimension(:), allocatable :: tmp
235 :
236 3 : n2 = 1_i4
237 :
238 3 : if (allocated(vec1)) then
239 1 : n1 = size(vec1)
240 : ! save vec1
241 1 : call move_alloc(vec1, tmp)
242 :
243 3 : allocate(vec1(n1+n2))
244 8 : vec1(1:n1) = tmp(1:n1)
245 1 : vec1(n1+1_i4) = sca2
246 : else
247 2 : n1 = 0_i4
248 :
249 2 : allocate(vec1(n2))
250 2 : vec1(1_i4) = sca2
251 : end if
252 :
253 3 : END SUBROUTINE append_i4_v_s
254 :
255 2 : SUBROUTINE append_i4_v_v(vec1, vec2)
256 :
257 : implicit none
258 :
259 : integer(i4), dimension(:), allocatable, intent(inout) :: vec1
260 : integer(i4), dimension(:), intent(in) :: vec2
261 :
262 : ! local variables
263 : integer(i4) :: n1, n2 ! length of vectors
264 2 : integer(i4), dimension(:), allocatable :: tmp
265 :
266 2 : n2 = size(vec2)
267 :
268 2 : if (allocated(vec1)) then
269 1 : n1 = size(vec1)
270 : ! save vec1
271 1 : call move_alloc(vec1, tmp)
272 :
273 3 : allocate(vec1(n1+n2))
274 6 : vec1(1:n1) = tmp(1:n1)
275 3 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
276 : else
277 1 : n1 = 0_i4
278 :
279 3 : allocate(vec1(n2))
280 6 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
281 : end if
282 :
283 5 : END SUBROUTINE append_i4_v_v
284 :
285 4 : SUBROUTINE append_i4_m_m(mat1, mat2, fill_value)
286 :
287 : implicit none
288 :
289 : integer(i4), dimension(:,:), allocatable, intent(inout) :: mat1
290 : integer(i4), dimension(:,:), intent(in) :: mat2
291 : integer(i4), optional, intent(in) :: fill_value
292 :
293 : ! local variables
294 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
295 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
296 4 : integer(i4), dimension(:,:), allocatable :: tmp
297 :
298 4 : m2 = size(mat2,1) ! rows
299 4 : n2 = size(mat2,2) ! columns
300 :
301 4 : if (allocated(mat1)) then
302 3 : m1 = size(mat1,1) ! rows
303 3 : n1 = size(mat1,2) ! columns
304 :
305 3 : if ((n1 /= n2) .and. .not. present(fill_value) ) then
306 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
307 0 : STOP
308 : end if
309 :
310 : ! save mat1
311 3 : call move_alloc(mat1, tmp)
312 :
313 3 : if ( n1 == n2 ) then
314 4 : allocate(mat1(m1+m2,n1))
315 7 : mat1(1:m1,:) = tmp(1:m1,:)
316 9 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
317 : end if
318 :
319 3 : if ( n1 > n2 ) then
320 4 : allocate(mat1(m1+m2,n1))
321 10 : mat1(1:m1,:) = tmp(1:m1,:)
322 7 : mat1(m1+1_i4:m1+m2, 1:n2) = mat2(1:m2,:)
323 4 : mat1(m1+1_i4:m1+m2,n2+1:n1) = fill_value
324 : end if
325 :
326 3 : if ( n1 < n2 ) then
327 4 : allocate(mat1(m1+m2,n2))
328 9 : mat1( 1:m1, 1:n1) = tmp(1:m1,:)
329 5 : mat1( 1:m1, n1+1:n2) = fill_value
330 10 : mat1(m1+1_i4:m1+m2, : ) = mat2(1:m2,:)
331 : end if
332 :
333 : else
334 1 : n1 = 0_i4
335 :
336 4 : allocate(mat1(m2,n2))
337 11 : mat1 = mat2
338 : end if
339 :
340 6 : END SUBROUTINE append_i4_m_m
341 :
342 2 : SUBROUTINE append_i4_3d(mat1, mat2, fill_value)
343 :
344 : implicit none
345 :
346 : integer(i4), dimension(:,:,:), allocatable, intent(inout) :: mat1
347 : integer(i4), dimension(:,:,:), intent(in) :: mat2
348 : integer(i4), optional, intent(in) :: fill_value
349 :
350 : ! local variables
351 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
352 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
353 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
354 2 : integer(i4), dimension(:,:,:), allocatable :: tmp
355 :
356 2 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_i4_3d'
357 :
358 2 : m2 = size(mat2,1) ! rows
359 2 : n2 = size(mat2,2) ! columns
360 2 : j2 = size(mat2,3) ! something else
361 :
362 2 : if (allocated(mat1)) then
363 1 : m1 = size(mat1,1) ! rows
364 1 : n1 = size(mat1,2) ! columns
365 1 : j1 = size(mat1,3) ! something else
366 :
367 1 : if ((n1 /= n2) .or. (j1 /= j2) ) then
368 : print*, 'append: size mismatch: dim 2 and 3 of matrix1 and matrix2 are unequal : ' &
369 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
370 0 : STOP
371 : end if
372 :
373 : ! save mat1
374 1 : call move_alloc(mat1, tmp)
375 :
376 5 : allocate(mat1(m1+m2,n1,j1))
377 15 : mat1(1:m1,:,:) = tmp(1:m1,:,:)
378 11 : mat1(m1+1_i4:m1+m2,:,:) = mat2(1:m2,:,:)
379 :
380 : else
381 :
382 5 : allocate(mat1(m2,n2,j2))
383 16 : mat1 = mat2
384 :
385 : end if
386 :
387 6 : END SUBROUTINE append_i4_3d
388 :
389 2 : SUBROUTINE append_i8_v_s(vec1, sca2)
390 :
391 : implicit none
392 :
393 : integer(i8), dimension(:), allocatable, intent(inout) :: vec1
394 : integer(i8), intent(in) :: sca2
395 :
396 : ! local variables
397 : integer(i4) :: n1, n2
398 2 : integer(i8), dimension(:), allocatable :: tmp
399 :
400 2 : n2 = 1_i4
401 :
402 2 : if (allocated(vec1)) then
403 1 : n1 = size(vec1)
404 : ! save vec1
405 1 : call move_alloc(vec1, tmp)
406 :
407 3 : allocate(vec1(n1+n2))
408 8 : vec1(1:n1) = tmp(1:n1)
409 1 : vec1(n1+1_i4) = sca2
410 : else
411 1 : n1 = 0_i4
412 :
413 1 : allocate(vec1(n2))
414 1 : vec1(1_i4) = sca2
415 : end if
416 :
417 4 : END SUBROUTINE append_i8_v_s
418 :
419 2 : SUBROUTINE append_i8_v_v(vec1, vec2)
420 :
421 : implicit none
422 :
423 : integer(i8), dimension(:), allocatable, intent(inout) :: vec1
424 : integer(i8), dimension(:), intent(in) :: vec2
425 :
426 : ! local variables
427 : integer(i4) :: n1, n2 ! length of vectors
428 2 : integer(i8), dimension(:), allocatable :: tmp
429 :
430 2 : n2 = size(vec2)
431 :
432 2 : if (allocated(vec1)) then
433 1 : n1 = size(vec1)
434 : ! save vec1
435 1 : call move_alloc(vec1, tmp)
436 :
437 3 : allocate(vec1(n1+n2))
438 6 : vec1(1:n1) = tmp(1:n1)
439 3 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
440 : else
441 1 : n1 = 0_i4
442 :
443 3 : allocate(vec1(n2))
444 6 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
445 : end if
446 :
447 4 : END SUBROUTINE append_i8_v_v
448 :
449 4 : SUBROUTINE append_i8_m_m(mat1, mat2, fill_value)
450 :
451 : implicit none
452 :
453 : integer(i8), dimension(:,:), allocatable, intent(inout) :: mat1
454 : integer(i8), dimension(:,:), intent(in) :: mat2
455 : integer(i8), optional, intent(in) :: fill_value
456 :
457 : ! local variables
458 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
459 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
460 4 : integer(i8), dimension(:,:), allocatable :: tmp
461 :
462 4 : m2 = size(mat2,1) ! rows
463 4 : n2 = size(mat2,2) ! columns
464 :
465 4 : if (allocated(mat1)) then
466 3 : m1 = size(mat1,1) ! rows
467 3 : n1 = size(mat1,2) ! columns
468 :
469 3 : if ((n1 /= n2) .and. .not. present(fill_value) ) then
470 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
471 0 : STOP
472 : end if
473 :
474 : ! save mat1
475 3 : call move_alloc(mat1, tmp)
476 :
477 3 : if ( n1 == n2 ) then
478 4 : allocate(mat1(m1+m2,n1))
479 7 : mat1(1:m1,:) = tmp(1:m1,:)
480 9 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
481 : end if
482 :
483 3 : if ( n1 > n2 ) then
484 4 : allocate(mat1(m1+m2,n1))
485 10 : mat1(1:m1,:) = tmp(1:m1,:)
486 7 : mat1(m1+1_i4:m1+m2, 1:n2) = mat2(1:m2,:)
487 4 : mat1(m1+1_i4:m1+m2,n2+1:n1) = fill_value
488 : end if
489 :
490 3 : if ( n1 < n2 ) then
491 4 : allocate(mat1(m1+m2,n2))
492 9 : mat1( 1:m1, 1:n1) = tmp(1:m1,:)
493 5 : mat1( 1:m1, n1+1:n2) = fill_value
494 10 : mat1(m1+1_i4:m1+m2, : ) = mat2(1:m2,:)
495 : end if
496 :
497 : else
498 1 : n1 = 0_i4
499 :
500 4 : allocate(mat1(m2,n2))
501 11 : mat1 = mat2
502 : end if
503 :
504 6 : END SUBROUTINE append_i8_m_m
505 :
506 2 : SUBROUTINE append_i8_3d(mat1, mat2, fill_value)
507 :
508 : implicit none
509 :
510 : integer(i8), dimension(:,:,:), allocatable, intent(inout) :: mat1
511 : integer(i8), dimension(:,:,:), intent(in) :: mat2
512 : integer(i8), optional, intent(in) :: fill_value
513 :
514 : ! local variables
515 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
516 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
517 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
518 2 : integer(i8), dimension(:,:,:), allocatable :: tmp
519 :
520 2 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_i8_3d'
521 :
522 2 : m2 = size(mat2,1) ! rows
523 2 : n2 = size(mat2,2) ! columns
524 2 : j2 = size(mat2,3) ! something else
525 :
526 2 : if (allocated(mat1)) then
527 1 : m1 = size(mat1,1) ! rows
528 1 : n1 = size(mat1,2) ! columns
529 1 : j1 = size(mat1,3) ! something else
530 :
531 1 : if ((n1 /= n2) .or. (j1 /= j2) ) then
532 : print*, 'append: size mismatch: dim 2 and 3 of matrix1 and matrix2 are unequal : ' &
533 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
534 0 : STOP
535 : end if
536 :
537 : ! save mat1
538 1 : call move_alloc(mat1, tmp)
539 :
540 5 : allocate(mat1(m1+m2,n1,j1))
541 15 : mat1(1:m1,:,:) = tmp(1:m1,:,:)
542 11 : mat1(m1+1_i4:m1+m2,:,:) = mat2(1:m2,:,:)
543 :
544 : else
545 :
546 5 : allocate(mat1(m2,n2,j2))
547 16 : mat1 = mat2
548 :
549 : end if
550 :
551 6 : END SUBROUTINE append_i8_3d
552 :
553 2 : SUBROUTINE append_sp_v_s(vec1, sca2)
554 :
555 : implicit none
556 :
557 : real(sp), dimension(:), allocatable, intent(inout) :: vec1
558 : real(sp), intent(in) :: sca2
559 :
560 : ! local variables
561 : integer(i4) :: n1, n2
562 2 : real(sp), dimension(:), allocatable :: tmp
563 :
564 2 : n2 = 1_i4
565 :
566 2 : if (allocated(vec1)) then
567 1 : n1 = size(vec1)
568 : ! save vec1
569 1 : call move_alloc(vec1, tmp)
570 :
571 3 : allocate(vec1(n1+n2))
572 8 : vec1(1:n1) = tmp(1:n1)
573 1 : vec1(n1+1_i4) = sca2
574 : else
575 1 : n1 = 0_i4
576 :
577 1 : allocate(vec1(n2))
578 1 : vec1(1_i4) = sca2
579 : end if
580 :
581 4 : END SUBROUTINE append_sp_v_s
582 :
583 2 : SUBROUTINE append_sp_v_v(vec1, vec2)
584 :
585 : implicit none
586 :
587 : real(sp), dimension(:), allocatable, intent(inout) :: vec1
588 : real(sp), dimension(:), intent(in) :: vec2
589 :
590 : ! local variables
591 : integer(i4) :: n1, n2 ! length of vectors
592 2 : real(sp), dimension(:), allocatable :: tmp
593 :
594 2 : n2 = size(vec2)
595 :
596 2 : if (allocated(vec1)) then
597 1 : n1 = size(vec1)
598 : ! save vec1
599 1 : call move_alloc(vec1, tmp)
600 :
601 3 : allocate(vec1(n1+n2))
602 6 : vec1(1:n1) = tmp(1:n1)
603 3 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
604 : else
605 1 : n1 = 0_i4
606 :
607 3 : allocate(vec1(n2))
608 6 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
609 : end if
610 :
611 4 : END SUBROUTINE append_sp_v_v
612 :
613 4 : SUBROUTINE append_sp_m_m(mat1, mat2, fill_value)
614 :
615 : implicit none
616 :
617 : real(sp), dimension(:,:), allocatable, intent(inout) :: mat1
618 : real(sp), dimension(:,:), intent(in) :: mat2
619 : real(sp), optional, intent(in) :: fill_value
620 :
621 : ! local variables
622 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
623 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
624 4 : real(sp), dimension(:,:), allocatable :: tmp
625 :
626 4 : m2 = size(mat2,1) ! rows
627 4 : n2 = size(mat2,2) ! columns
628 :
629 4 : if (allocated(mat1)) then
630 3 : m1 = size(mat1,1) ! rows
631 3 : n1 = size(mat1,2) ! columns
632 :
633 3 : if ((n1 /= n2) .and. .not. present(fill_value) ) then
634 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
635 0 : STOP
636 : end if
637 :
638 : ! save mat1
639 3 : call move_alloc(mat1, tmp)
640 :
641 3 : if ( n1 == n2 ) then
642 4 : allocate(mat1(m1+m2,n1))
643 7 : mat1(1:m1,:) = tmp(1:m1,:)
644 9 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
645 : end if
646 :
647 3 : if ( n1 > n2 ) then
648 4 : allocate(mat1(m1+m2,n1))
649 10 : mat1(1:m1,:) = tmp(1:m1,:)
650 7 : mat1(m1+1_i4:m1+m2, 1:n2) = mat2(1:m2,:)
651 4 : mat1(m1+1_i4:m1+m2,n2+1:n1) = fill_value
652 : end if
653 :
654 3 : if ( n1 < n2 ) then
655 4 : allocate(mat1(m1+m2,n2))
656 9 : mat1( 1:m1, 1:n1) = tmp(1:m1,:)
657 5 : mat1( 1:m1, n1+1:n2) = fill_value
658 10 : mat1(m1+1_i4:m1+m2, : ) = mat2(1:m2,:)
659 : end if
660 :
661 : else
662 1 : n1 = 0_i4
663 :
664 4 : allocate(mat1(m2,n2))
665 11 : mat1 = mat2
666 : end if
667 :
668 6 : END SUBROUTINE append_sp_m_m
669 :
670 2 : SUBROUTINE append_sp_3d(mat1, mat2, fill_value)
671 :
672 : implicit none
673 :
674 : real(sp), dimension(:,:,:), allocatable, intent(inout) :: mat1
675 : real(sp), dimension(:,:,:), intent(in) :: mat2
676 : real(sp), optional, intent(in) :: fill_value
677 :
678 : ! local variables
679 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
680 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
681 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
682 2 : real(sp), dimension(:,:,:), allocatable :: tmp
683 :
684 2 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_sp_3d'
685 :
686 2 : m2 = size(mat2,1) ! rows
687 2 : n2 = size(mat2,2) ! columns
688 2 : j2 = size(mat2,3) ! something else
689 :
690 2 : if (allocated(mat1)) then
691 1 : m1 = size(mat1,1) ! rows
692 1 : n1 = size(mat1,2) ! columns
693 1 : j1 = size(mat1,3) ! something else
694 :
695 1 : if ((n1 /= n2) .or. (j1 /= j2) ) then
696 : print*, 'append: size mismatch: dim 2 and 3 of matrix1 and matrix2 are unequal : ' &
697 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
698 0 : STOP
699 : end if
700 :
701 : ! save mat1
702 1 : call move_alloc(mat1, tmp)
703 :
704 5 : allocate(mat1(m1+m2,n1,j1))
705 15 : mat1(1:m1,:,:) = tmp(1:m1,:,:)
706 11 : mat1(m1+1_i4:m1+m2,:,:) = mat2(1:m2,:,:)
707 :
708 : else
709 :
710 5 : allocate(mat1(m2,n2,j2))
711 16 : mat1 = mat2
712 :
713 : end if
714 :
715 6 : END SUBROUTINE append_sp_3d
716 :
717 52 : SUBROUTINE append_dp_v_s(vec1, sca2)
718 :
719 : implicit none
720 :
721 : real(dp), dimension(:), allocatable, intent(inout) :: vec1
722 : real(dp), intent(in) :: sca2
723 :
724 : ! local variables
725 : integer(i4) :: n1, n2
726 52 : real(dp), dimension(:), allocatable :: tmp
727 :
728 52 : n2 = 1_i4
729 :
730 52 : if (allocated(vec1)) then
731 49 : n1 = size(vec1)
732 : ! save vec1
733 49 : call move_alloc(vec1, tmp)
734 :
735 147 : allocate(vec1(n1+n2))
736 429 : vec1(1:n1) = tmp(1:n1)
737 49 : vec1(n1+1_i4) = sca2
738 : else
739 3 : n1 = 0_i4
740 :
741 3 : allocate(vec1(n2))
742 3 : vec1(1_i4) = sca2
743 : end if
744 :
745 54 : END SUBROUTINE append_dp_v_s
746 :
747 2 : SUBROUTINE append_dp_v_v(vec1, vec2)
748 :
749 : implicit none
750 :
751 : real(dp), dimension(:), allocatable, intent(inout) :: vec1
752 : real(dp), dimension(:), intent(in) :: vec2
753 :
754 : ! local variables
755 : integer(i4) :: n1, n2 ! length of vectors
756 2 : real(dp), dimension(:), allocatable :: tmp
757 :
758 2 : n2 = size(vec2)
759 :
760 2 : if (allocated(vec1)) then
761 1 : n1 = size(vec1)
762 : ! save vec1
763 1 : call move_alloc(vec1, tmp)
764 :
765 3 : allocate(vec1(n1+n2))
766 6 : vec1(1:n1) = tmp(1:n1)
767 3 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
768 : else
769 1 : n1 = 0_i4
770 :
771 3 : allocate(vec1(n2))
772 6 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
773 : end if
774 :
775 54 : END SUBROUTINE append_dp_v_v
776 :
777 32 : SUBROUTINE append_dp_m_m(mat1, mat2, fill_value, idim)
778 :
779 : implicit none
780 :
781 : real(dp), dimension(:,:), allocatable, intent(inout) :: mat1
782 : real(dp), dimension(:,:), intent(in) :: mat2
783 : real(dp), optional, intent(in) :: fill_value
784 : integer(i4), optional, intent(in) :: idim
785 :
786 : ! local variables
787 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
788 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
789 : integer(i4) :: dd
790 32 : real(dp), dimension(:,:), allocatable :: tmp
791 :
792 32 : dd = 1
793 32 : if (present(idim)) dd = idim
794 32 : if (dd > 2) then
795 0 : print*, 'append: dd is : (',dd,') and greater than number of dimensions : 2'
796 0 : STOP
797 : end if
798 :
799 :
800 32 : m2 = size(mat2,1) ! rows
801 32 : n2 = size(mat2,2) ! columns
802 :
803 32 : if (allocated(mat1)) then
804 29 : m1 = size(mat1,1) ! rows
805 29 : n1 = size(mat1,2) ! columns
806 :
807 29 : if (present(idim)) then
808 2 : if (dd == 1) then
809 1 : if (n1 /= n2) then
810 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
811 0 : STOP
812 : end if
813 :
814 : ! save mat1
815 1 : call move_alloc(mat1, tmp)
816 :
817 4 : allocate(mat1(m1+m2,n1))
818 16 : mat1(1:m1,:) = tmp(1:m1,:)
819 10 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
820 :
821 1 : else if (dd == 2) then
822 :
823 1 : if (m1 /= m2) then
824 0 : print*, 'append: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
825 0 : STOP
826 : end if
827 :
828 : ! save mat1
829 1 : call move_alloc(mat1, tmp)
830 :
831 4 : allocate(mat1(m1,n1 + n2))
832 10 : mat1(:,1:n1) = tmp(:,1:n1)
833 7 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
834 :
835 : end if
836 :
837 : else
838 :
839 27 : if ((n1 /= n2) .and. .not. present(fill_value) ) then
840 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
841 0 : STOP
842 : end if
843 :
844 : ! save mat1
845 27 : call move_alloc(mat1, tmp)
846 :
847 27 : if ( n1 == n2 ) then
848 100 : allocate(mat1(m1+m2,n1))
849 540436 : mat1(1:m1,:) = tmp(1:m1,:)
850 84186 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
851 : end if
852 :
853 27 : if ( n1 > n2 ) then
854 4 : allocate(mat1(m1+m2,n1))
855 10 : mat1(1:m1,:) = tmp(1:m1,:)
856 7 : mat1(m1+1_i4:m1+m2, 1:n2) = mat2(1:m2,:)
857 4 : mat1(m1+1_i4:m1+m2,n2+1:n1) = fill_value
858 : end if
859 :
860 27 : if ( n1 < n2 ) then
861 4 : allocate(mat1(m1+m2,n2))
862 9 : mat1( 1:m1, 1:n1) = tmp(1:m1,:)
863 5 : mat1( 1:m1, n1+1:n2) = fill_value
864 10 : mat1(m1+1_i4:m1+m2, : ) = mat2(1:m2,:)
865 : end if
866 : end if
867 :
868 : else
869 12 : allocate(mat1(m2,n2))
870 1803 : mat1 = mat2
871 : end if
872 :
873 34 : END SUBROUTINE append_dp_m_m
874 :
875 5 : SUBROUTINE append_dp_3d(mat1, mat2, fill_value, idim)
876 :
877 : implicit none
878 :
879 : real(dp), dimension(:,:,:), allocatable, intent(inout) :: mat1
880 : real(dp), dimension(:,:,:), intent(in) :: mat2
881 : real(dp), optional, intent(in) :: fill_value
882 : integer(i4), optional, intent(in) :: idim
883 :
884 : ! local variables
885 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
886 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
887 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
888 : integer(i4) :: dd ! dimension it is appended along
889 5 : real(dp), dimension(:,:,:), allocatable :: tmp
890 :
891 5 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_dp_3d'
892 :
893 5 : dd = 1_i4
894 5 : if (present(idim)) dd = idim
895 5 : if (dd > 3) then
896 0 : print*, 'append: dd is : (',dd,') and greater than number of dimensions : 3'
897 0 : STOP
898 : end if
899 :
900 5 : m2 = size(mat2,1) ! rows
901 5 : n2 = size(mat2,2) ! columns
902 5 : j2 = size(mat2,3) ! something else
903 :
904 5 : if (allocated(mat1)) then
905 4 : m1 = size(mat1,1) ! rows
906 4 : n1 = size(mat1,2) ! columns
907 4 : j1 = size(mat1,3) ! something else
908 :
909 4 : if (dd == 1) then
910 2 : if ((n1 /= n2) .or. (j1 /= j2) ) then
911 : print*, 'append: size mismatch: dim 2 and 3 of matrix1 and matrix2 are unequal : ' &
912 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
913 0 : STOP
914 : end if
915 :
916 : ! save mat1
917 2 : call move_alloc(mat1, tmp)
918 :
919 10 : allocate(mat1(m1+m2,n1,j1))
920 22 : mat1(1:m1,:,:) = tmp(1:m1,:,:)
921 22 : mat1(m1+1_i4:m1+m2,:,:) = mat2(1:m2,:,:)
922 2 : else if (dd == 2) then
923 1 : if ((m1 /= m2) .or. (j1 /= j2) ) then
924 : print*, 'append: size mismatch: dim 1 and 3 of matrix1 and matrix2 are unequal : ' &
925 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
926 0 : STOP
927 : end if
928 :
929 : ! save mat1
930 1 : call move_alloc(mat1, tmp)
931 :
932 5 : allocate(mat1(m1,n1 + n2,j1))
933 15 : mat1(:,1:n1,:) = tmp(:,1:n1,:)
934 9 : mat1(:,n1+1_i4:n1+n2,:) = mat2(:,1:n2,:)
935 1 : else if (dd == 3) then
936 1 : if ((m1 /= m2) .or. (n1 /= n2) ) then
937 : print*, 'append: size mismatch: dim 1 and 2 of matrix1 and matrix2 are unequal : ' &
938 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
939 0 : STOP
940 : end if
941 :
942 : ! save mat1
943 1 : call move_alloc(mat1, tmp)
944 :
945 5 : allocate(mat1(m1,n1,j1 + j2))
946 15 : mat1(:,:,1:j1) = tmp(:,:,1:j1)
947 8 : mat1(:,:,j1+1_i4:j1+j2) = mat2(:,:,1:j2)
948 : end if
949 :
950 : else
951 :
952 5 : allocate(mat1(m2,n2,j2))
953 12 : mat1 = mat2
954 :
955 : end if
956 :
957 37 : END SUBROUTINE append_dp_3d
958 :
959 6 : SUBROUTINE append_dp_4d(mat1, mat2, fill_value, idim)
960 :
961 : implicit none
962 :
963 : real(dp), dimension(:,:,:,:), allocatable, intent(inout) :: mat1
964 : real(dp), dimension(:,:,:,:), intent(in) :: mat2
965 : real(dp), optional, intent(in) :: fill_value
966 : integer(i4), optional, intent(in) :: idim
967 :
968 : ! local variables
969 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
970 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
971 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
972 : integer(i4) :: i1, i2 ! dim4 of matrixes: something else
973 : integer(i4) :: dd ! dimension it is appended along
974 6 : real(dp), allocatable :: tmp(:,:,:,:)
975 :
976 6 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_dp_3d'
977 :
978 6 : dd = 1_i4
979 6 : if (present(idim)) dd = idim
980 6 : if (dd > 4) then
981 0 : print*, 'append: dd is : (',dd,') and greater than number of dimensions : 4'
982 0 : STOP
983 : end if
984 :
985 6 : m2 = size(mat2,1) ! rows
986 6 : n2 = size(mat2,2) ! columns
987 6 : j2 = size(mat2,3) ! something else
988 6 : i2 = size(mat2,4) ! something else
989 :
990 6 : if (allocated(mat1)) then
991 5 : m1 = size(mat1,1) ! rows
992 5 : n1 = size(mat1,2) ! columns
993 5 : j1 = size(mat1,3) ! something else
994 5 : i1 = size(mat1,4) ! something else
995 :
996 5 : if (dd == 1) then
997 2 : if ((n1 /= n2) .or. (j1 /= j2) .or. (i1 /= i2)) then
998 : print*, 'append: size mismatch: dim 2, 3, and 4 of matrix1 and matrix2 are unequal : ' &
999 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1000 0 : STOP
1001 : end if
1002 :
1003 : ! save mat1
1004 2 : call move_alloc(mat1, tmp)
1005 :
1006 12 : allocate(mat1(m1+m2,n1,j1,i1))
1007 24 : mat1(1:m1,:,:,:) = tmp(1:m1,:,:,:)
1008 24 : mat1(m1+1_i4:m1+m2,:,:,:) = mat2(1:m2,:,:,:)
1009 3 : else if (dd == 2) then
1010 1 : if ((m1 /= m2) .or. (j1 /= j2) .or. (i1 /= i2)) then
1011 : print*, 'append: size mismatch: dim 1, 3, and 4 of matrix1 and matrix2 are unequal : ' &
1012 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1013 0 : STOP
1014 : end if
1015 :
1016 : ! save mat1
1017 1 : call move_alloc(mat1, tmp)
1018 :
1019 6 : allocate(mat1(m1,n1 + n2,j1,i1))
1020 16 : mat1(:,1:n1,:,:) = tmp(:,1:n1,:,:)
1021 10 : mat1(:,n1+1_i4:n1+n2,:,:) = mat2(:,1:n2,:,:)
1022 2 : else if (dd == 3) then
1023 1 : if ((m1 /= m2) .or. (n1 /= n2) .or. (i1 /= i2)) then
1024 : print*, 'append: size mismatch: dim 1, 2, and 4 of matrix1 and matrix2 are unequal : ' &
1025 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1026 0 : STOP
1027 : end if
1028 :
1029 : ! save mat1
1030 6 : allocate(tmp(m1,n1,j1,i1))
1031 17 : tmp=mat1
1032 1 : deallocate(mat1)
1033 :
1034 6 : allocate(mat1(m1,n1,j1 + j2,i1))
1035 16 : mat1(:,:,1:j1,:) = tmp(:,:,1:j1,:)
1036 9 : mat1(:,:,j1+1_i4:j1+j2,:) = mat2(:,:,1:j2,:)
1037 1 : else if (dd == 4) then
1038 1 : if ((m1 /= m2) .or. (n1 /= n2) .or. (j1 /= j2)) then
1039 : print*, 'append: size mismatch: dim 1, 2, and 3 of matrix1 and matrix2 are unequal : ' &
1040 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1041 0 : STOP
1042 : end if
1043 :
1044 : ! save mat1
1045 1 : call move_alloc(mat1, tmp)
1046 :
1047 6 : allocate(mat1(m1,n1,j1,i1 + i2))
1048 23 : mat1(:,:,:,1:i1) = tmp(:,:,:,1:i1)
1049 23 : mat1(:,:,:,i1+1_i4:i1+i2) = mat2(:,:,:,1:i2)
1050 : end if
1051 :
1052 : else
1053 :
1054 6 : allocate(mat1(m2,n2,j2,i2))
1055 13 : mat1 = mat2
1056 :
1057 : end if
1058 :
1059 11 : END SUBROUTINE append_dp_4d
1060 :
1061 220 : SUBROUTINE append_char_v_s(vec1, sca2)
1062 :
1063 : implicit none
1064 :
1065 : character(len=*), dimension(:), allocatable, intent(inout) :: vec1
1066 : character(len=*), intent(in) :: sca2
1067 :
1068 : ! local variables
1069 : integer(i4) :: n1, n2
1070 220 : character(len(vec1)), dimension(:), allocatable :: tmp
1071 :
1072 220 : n2 = 1_i4
1073 :
1074 220 : if (allocated(vec1)) then
1075 216 : n1 = size(vec1)
1076 : ! save mat1
1077 216 : call move_alloc(vec1, tmp)
1078 :
1079 648 : allocate(vec1(n1+n2))
1080 1525 : vec1(1:n1) = tmp(1:n1)
1081 216 : vec1(n1+1_i4) = sca2
1082 : else
1083 4 : n1 = 0_i4
1084 :
1085 4 : allocate(vec1(n2))
1086 4 : vec1(1_i4) = sca2
1087 : end if
1088 :
1089 226 : END SUBROUTINE append_char_v_s
1090 :
1091 2 : SUBROUTINE append_char_v_v(vec1, vec2)
1092 :
1093 : character(len=*), dimension(:), allocatable, intent(inout) :: vec1
1094 : character(len=*), dimension(:), intent(in) :: vec2
1095 :
1096 : ! local variables
1097 : integer(i4) :: n1, n2
1098 2 : character(len(vec1)), dimension(:), allocatable :: tmp
1099 :
1100 2 : n2 = size(vec2)
1101 :
1102 2 : if (allocated(vec1)) then
1103 1 : n1 = size(vec1)
1104 : ! save mat1
1105 1 : call move_alloc(vec1, tmp)
1106 :
1107 3 : allocate(vec1(n1+n2))
1108 6 : vec1(1:n1) = tmp(1:n1)
1109 3 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
1110 : else
1111 1 : n1 = 0_i4
1112 :
1113 3 : allocate(vec1(n2))
1114 6 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
1115 : end if
1116 :
1117 222 : END SUBROUTINE append_char_v_v
1118 :
1119 4 : SUBROUTINE append_char_m_m(mat1, mat2, fill_value)
1120 :
1121 : implicit none
1122 :
1123 : character(len=*), dimension(:,:), allocatable, intent(inout) :: mat1
1124 : character(len=*), dimension(:,:), intent(in) :: mat2
1125 : character(len=*), optional, intent(in) :: fill_value
1126 :
1127 : ! local variables
1128 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
1129 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
1130 4 : character(len(mat1)), dimension(:,:), allocatable :: tmp
1131 :
1132 4 : m2 = size(mat2,1) ! rows
1133 4 : n2 = size(mat2,2) ! columns
1134 :
1135 4 : if (allocated(mat1)) then
1136 3 : m1 = size(mat1,1) ! rows
1137 3 : n1 = size(mat1,2) ! columns
1138 :
1139 3 : if ((n1 /= n2) .and. .not. present(fill_value) ) then
1140 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1141 0 : STOP
1142 : end if
1143 :
1144 : ! save mat1
1145 3 : call move_alloc(mat1, tmp)
1146 :
1147 3 : if ( n1 == n2 ) then
1148 4 : allocate(mat1(m1+m2,n1))
1149 7 : mat1(1:m1,:) = tmp(1:m1,:)
1150 9 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
1151 : end if
1152 :
1153 3 : if ( n1 > n2 ) then
1154 4 : allocate(mat1(m1+m2,n1))
1155 13 : mat1(1:m1,:) = tmp(1:m1,:)
1156 5 : mat1(m1+1_i4:m1+m2, 1:n2) = mat2(1:m2,:)
1157 6 : mat1(m1+1_i4:m1+m2,n2+1:n1) = fill_value
1158 : end if
1159 :
1160 3 : if ( n1 < n2 ) then
1161 4 : allocate(mat1(m1+m2,n2))
1162 19 : mat1( 1:m1, 1:n1) = tmp(1:m1,:)
1163 11 : mat1( 1:m1, n1+1:n2) = fill_value
1164 10 : mat1(m1+1_i4:m1+m2, : ) = mat2(1:m2,:)
1165 : end if
1166 :
1167 : else
1168 1 : n1 = 0_i4
1169 :
1170 4 : allocate(mat1(m2,n2))
1171 8 : mat1 = mat2
1172 : end if
1173 :
1174 6 : END SUBROUTINE append_char_m_m
1175 :
1176 2 : SUBROUTINE append_char_3d(mat1, mat2, fill_value)
1177 :
1178 : implicit none
1179 :
1180 : character(len=*), dimension(:,:,:), allocatable, intent(inout) :: mat1
1181 : character(len=*), dimension(:,:,:), intent(in) :: mat2
1182 : character(len=*), optional, intent(in) :: fill_value
1183 :
1184 : ! local variables
1185 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
1186 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
1187 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
1188 2 : character(len(mat1)), dimension(:,:,:), allocatable :: tmp
1189 :
1190 2 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_i8_3d'
1191 :
1192 2 : m2 = size(mat2,1) ! rows
1193 2 : n2 = size(mat2,2) ! columns
1194 2 : j2 = size(mat2,3) ! something else
1195 :
1196 2 : if (allocated(mat1)) then
1197 1 : m1 = size(mat1,1) ! rows
1198 1 : n1 = size(mat1,2) ! columns
1199 1 : j1 = size(mat1,3) ! something else
1200 :
1201 1 : if ((n1 /= n2) .or. (j1 /= j2) ) then
1202 : print*, 'append: size mismatch: dim 2 and 3 of matrix1 and matrix2 are unequal : ' &
1203 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
1204 0 : STOP
1205 : end if
1206 :
1207 : ! save mat1
1208 1 : call move_alloc(mat1, tmp)
1209 :
1210 5 : allocate(mat1(m1+m2,n1,j1))
1211 15 : mat1(1:m1,:,:) = tmp(1:m1,:,:)
1212 11 : mat1(m1+1_i4:m1+m2,:,:) = mat2(1:m2,:,:)
1213 :
1214 : else
1215 :
1216 5 : allocate(mat1(m2,n2,j2))
1217 16 : mat1 = mat2
1218 :
1219 : end if
1220 :
1221 6 : END SUBROUTINE append_char_3d
1222 :
1223 2 : SUBROUTINE append_lgt_v_s(vec1, sca2)
1224 :
1225 : implicit none
1226 :
1227 : logical, dimension(:), allocatable, intent(inout) :: vec1
1228 : logical, intent(in) :: sca2
1229 :
1230 : ! local variables
1231 : integer(i4) :: n1, n2
1232 2 : logical, dimension(:), allocatable :: tmp
1233 :
1234 2 : n2 = 1_i4
1235 :
1236 2 : if (allocated(vec1)) then
1237 1 : n1 = size(vec1)
1238 : ! save mat1
1239 1 : call move_alloc(vec1, tmp)
1240 :
1241 3 : allocate(vec1(n1+n2))
1242 8 : vec1(1:n1) = tmp(1:n1)
1243 1 : vec1(n1+1_i4) = sca2
1244 : else
1245 1 : n1 = 0_i4
1246 :
1247 1 : allocate(vec1(n2))
1248 1 : vec1(1_i4) = sca2
1249 : end if
1250 :
1251 4 : END SUBROUTINE append_lgt_v_s
1252 :
1253 2 : SUBROUTINE append_lgt_v_v(vec1, vec2)
1254 :
1255 : implicit none
1256 :
1257 : logical, dimension(:), allocatable, intent(inout) :: vec1
1258 : logical, dimension(:), intent(in) :: vec2
1259 :
1260 : ! local variables
1261 : integer(i4) :: n1, n2 ! length of vectors
1262 2 : logical, dimension(:), allocatable :: tmp
1263 :
1264 2 : n2 = size(vec2)
1265 :
1266 2 : if (allocated(vec1)) then
1267 1 : n1 = size(vec1)
1268 : ! save mat1
1269 1 : call move_alloc(vec1, tmp)
1270 :
1271 3 : allocate(vec1(n1+n2))
1272 6 : vec1(1:n1) = tmp(1:n1)
1273 3 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
1274 : else
1275 1 : n1 = 0_i4
1276 :
1277 3 : allocate(vec1(n2))
1278 6 : vec1(n1+1_i4:n1+n2) = vec2(1:n2)
1279 : end if
1280 :
1281 4 : END SUBROUTINE append_lgt_v_v
1282 :
1283 6 : SUBROUTINE append_lgt_m_m(mat1, mat2, fill_value, idim)
1284 :
1285 : implicit none
1286 :
1287 : logical, dimension(:,:), allocatable, intent(inout) :: mat1
1288 : logical, dimension(:,:), intent(in) :: mat2
1289 : logical, optional, intent(in) :: fill_value
1290 : integer(i4), optional, intent(in) :: idim
1291 :
1292 : ! local variables
1293 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
1294 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
1295 : integer(i4) :: dd
1296 6 : logical, dimension(:,:), allocatable :: tmp
1297 :
1298 6 : dd = 1
1299 6 : if (present(idim)) dd = idim
1300 6 : if (dd > 3) then
1301 0 : print*, 'append: dd is : (',dd,') and greater than number of dimensions : 3'
1302 0 : STOP
1303 : end if
1304 :
1305 6 : m2 = size(mat2,1) ! rows
1306 6 : n2 = size(mat2,2) ! columns
1307 :
1308 6 : if (allocated(mat1)) then
1309 5 : m1 = size(mat1,1) ! rows
1310 5 : n1 = size(mat1,2) ! columns
1311 :
1312 5 : if (present(idim)) then
1313 2 : if (dd == 1) then
1314 1 : if (n1 /= n2) then
1315 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1316 0 : STOP
1317 : end if
1318 :
1319 : ! save mat1
1320 1 : call move_alloc(mat1, tmp)
1321 :
1322 4 : allocate(mat1(m1+m2,n1))
1323 16 : mat1(1:m1,:) = tmp(1:m1,:)
1324 10 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
1325 :
1326 1 : else if (dd == 2) then
1327 :
1328 1 : if (m1 /= m2) then
1329 0 : print*, 'append: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1330 0 : STOP
1331 : end if
1332 :
1333 : ! save mat1
1334 1 : call move_alloc(mat1, tmp)
1335 :
1336 4 : allocate(mat1(m1,n1 + n2))
1337 10 : mat1(:,1:n1) = tmp(:,1:n1)
1338 7 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1339 :
1340 : end if
1341 :
1342 : else
1343 3 : if ( (n1 /= n2) .and. .not. present(fill_value) ) then
1344 0 : print*, 'append: columns of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1345 0 : STOP
1346 : end if
1347 :
1348 : ! save mat1
1349 3 : call move_alloc(mat1, tmp)
1350 :
1351 3 : if ( n1 == n2 ) then
1352 4 : allocate(mat1(m1+m2,n1))
1353 7 : mat1(1:m1,:) = tmp(1:m1,:)
1354 9 : mat1(m1+1_i4:m1+m2,:) = mat2(1:m2,:)
1355 : end if
1356 :
1357 3 : if ( n1 > n2 ) then
1358 4 : allocate(mat1(m1+m2,n1))
1359 10 : mat1(1:m1,:) = tmp(1:m1,:)
1360 7 : mat1(m1+1_i4:m1+m2, 1:n2) = mat2(1:m2,:)
1361 4 : mat1(m1+1_i4:m1+m2,n2+1:n1) = fill_value
1362 : end if
1363 :
1364 3 : if ( n1 < n2 ) then
1365 4 : allocate(mat1(m1+m2,n2))
1366 9 : mat1( 1:m1, 1:n1) = tmp(1:m1,:)
1367 5 : mat1( 1:m1, n1+1:n2) = fill_value
1368 10 : mat1(m1+1_i4:m1+m2, : ) = mat2(1:m2,:)
1369 : end if
1370 : end if
1371 :
1372 : else
1373 1 : n1 = 0_i4
1374 :
1375 4 : allocate(mat1(m2,n2))
1376 11 : mat1 = mat2
1377 : end if
1378 :
1379 8 : END SUBROUTINE append_lgt_m_m
1380 :
1381 5 : SUBROUTINE append_lgt_3d(mat1, mat2, fill_value, idim)
1382 :
1383 : implicit none
1384 :
1385 : logical, dimension(:,:,:), allocatable, intent(inout) :: mat1
1386 : logical, dimension(:,:,:), intent(in) :: mat2
1387 : logical, optional, intent(in) :: fill_value
1388 : integer(i4), optional :: idim
1389 :
1390 : ! local variables
1391 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
1392 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
1393 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
1394 : integer(i4) :: dd ! index of dimension to append
1395 5 : logical, dimension(:,:,:), allocatable :: tmp
1396 :
1397 5 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_i8_3d'
1398 :
1399 5 : dd = 1
1400 5 : if (present(idim)) dd = idim
1401 5 : if (dd > 3) then
1402 0 : print*, 'append: dd is : (',dd,') and greater than number of dimensions : 3'
1403 0 : STOP
1404 : end if
1405 :
1406 5 : m2 = size(mat2,1) ! rows
1407 5 : n2 = size(mat2,2) ! columns
1408 5 : j2 = size(mat2,3) ! something else
1409 :
1410 5 : if (allocated(mat1)) then
1411 4 : m1 = size(mat1,1) ! rows
1412 4 : n1 = size(mat1,2) ! columns
1413 4 : j1 = size(mat1,3) ! something else
1414 :
1415 4 : if (dd == 1) then
1416 2 : if ((n1 /= n2) .or. (j1 /= j2) ) then
1417 : print*, 'append: size mismatch: dim 2 and 3 of matrix1 and matrix2 are unequal : ' &
1418 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
1419 0 : STOP
1420 : end if
1421 :
1422 : ! save mat1
1423 2 : call move_alloc(mat1, tmp)
1424 :
1425 10 : allocate(mat1(m1+m2,n1,j1))
1426 22 : mat1(1:m1,:,:) = tmp(1:m1,:,:)
1427 22 : mat1(m1+1_i4:m1+m2,:,:) = mat2(1:m2,:,:)
1428 2 : else if (dd == 2) then
1429 1 : if ((m1 /= m2) .or. (j1 /= j2) ) then
1430 : print*, 'append: size mismatch: dim 1 and 3 of matrix1 and matrix2 are unequal : ' &
1431 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
1432 0 : STOP
1433 : end if
1434 :
1435 : ! save mat1
1436 1 : call move_alloc(mat1, tmp)
1437 :
1438 5 : allocate(mat1(m1,n1 + n2,j1))
1439 15 : mat1(:,1:n1,:) = tmp(:,1:n1,:)
1440 9 : mat1(:,n1+1_i4:n1+n2,:) = mat2(:,1:n2,:)
1441 1 : else if (dd == 3) then
1442 1 : if ((m1 /= m2) .or. (n1 /= n2) ) then
1443 : print*, 'append: size mismatch: dim 1 and 2 of matrix1 and matrix2 are unequal : ' &
1444 0 : // '(',m1,',',n1,',',j1,') and (',m2,',',n2,',',j2,')'
1445 0 : STOP
1446 : end if
1447 :
1448 : ! save mat1
1449 1 : call move_alloc(mat1, tmp)
1450 :
1451 5 : allocate(mat1(m1,n1,j1 + j2))
1452 15 : mat1(:,:,1:j1) = tmp(:,:,1:j1)
1453 8 : mat1(:,:,j1+1_i4:j1+j2) = mat2(:,:,1:j2)
1454 : end if
1455 :
1456 : else
1457 :
1458 5 : allocate(mat1(m2,n2,j2))
1459 12 : mat1 = mat2
1460 :
1461 : end if
1462 :
1463 11 : END SUBROUTINE append_lgt_3d
1464 :
1465 6 : SUBROUTINE append_lgt_4d(mat1, mat2, fill_value, idim)
1466 :
1467 : implicit none
1468 :
1469 : logical, dimension(:,:,:,:), allocatable, intent(inout) :: mat1
1470 : logical, dimension(:,:,:,:), intent(in) :: mat2
1471 : logical, optional, intent(in) :: fill_value
1472 : integer(i4), optional, intent(in) :: idim
1473 :
1474 : ! local variables
1475 : integer(i4) :: m1, m2 ! dim1 of matrixes: rows
1476 : integer(i4) :: n1, n2 ! dim2 of matrixes: columns
1477 : integer(i4) :: j1, j2 ! dim3 of matrixes: something else
1478 : integer(i4) :: i1, i2 ! dim4 of matrixes: something else
1479 : integer(i4) :: dd ! dimension it is appended along
1480 6 : logical, allocatable :: tmp(:,:,:,:)
1481 :
1482 6 : if (present(fill_value)) print*, '***warning: fill_value is ignored in append_lgt_4d'
1483 :
1484 6 : dd = 1_i4
1485 6 : if (present(idim)) dd = idim
1486 6 : if (dd > 4) then
1487 0 : print*, 'append: dd is : (',dd,') and greater than number of dimensions : 4'
1488 0 : STOP
1489 : end if
1490 :
1491 6 : m2 = size(mat2,1) ! rows
1492 6 : n2 = size(mat2,2) ! columns
1493 6 : j2 = size(mat2,3) ! something else
1494 6 : i2 = size(mat2,4) ! something else
1495 :
1496 6 : if (allocated(mat1)) then
1497 5 : m1 = size(mat1,1) ! rows
1498 5 : n1 = size(mat1,2) ! columns
1499 5 : j1 = size(mat1,3) ! something else
1500 5 : i1 = size(mat1,4) ! something else
1501 :
1502 5 : if (dd == 1) then
1503 2 : if ((n1 /= n2) .or. (j1 /= j2) .or. (i1 /= i2)) then
1504 : print*, 'append: size mismatch: dim 2, 3, and 4 of matrix1 and matrix2 are unequal : ' &
1505 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1506 0 : STOP
1507 : end if
1508 :
1509 : ! save mat1
1510 2 : call move_alloc(mat1, tmp)
1511 :
1512 12 : allocate(mat1(m1+m2,n1,j1,i1))
1513 24 : mat1(1:m1,:,:,:) = tmp(1:m1,:,:,:)
1514 24 : mat1(m1+1_i4:m1+m2,:,:,:) = mat2(1:m2,:,:,:)
1515 3 : else if (dd == 2) then
1516 1 : if ((m1 /= m2) .or. (j1 /= j2) .or. (i1 /= i2)) then
1517 : print*, 'append: size mismatch: dim 1, 3, and 4 of matrix1 and matrix2 are unequal : ' &
1518 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1519 0 : STOP
1520 : end if
1521 :
1522 : ! save mat1
1523 1 : call move_alloc(mat1, tmp)
1524 :
1525 6 : allocate(mat1(m1,n1 + n2,j1,i1))
1526 16 : mat1(:,1:n1,:,:) = tmp(:,1:n1,:,:)
1527 10 : mat1(:,n1+1_i4:n1+n2,:,:) = mat2(:,1:n2,:,:)
1528 2 : else if (dd == 3) then
1529 1 : if ((m1 /= m2) .or. (n1 /= n2) .or. (i1 /= i2)) then
1530 : print*, 'append: size mismatch: dim 1, 2, and 4 of matrix1 and matrix2 are unequal : ' &
1531 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1532 0 : STOP
1533 : end if
1534 :
1535 : ! save mat1
1536 1 : call move_alloc(mat1, tmp)
1537 :
1538 6 : allocate(mat1(m1,n1,j1 + j2,i1))
1539 16 : mat1(:,:,1:j1,:) = tmp(:,:,1:j1,:)
1540 9 : mat1(:,:,j1+1_i4:j1+j2,:) = mat2(:,:,1:j2,:)
1541 1 : else if (dd == 4) then
1542 1 : if ((m1 /= m2) .or. (n1 /= n2) .or. (j1 /= j2)) then
1543 : print*, 'append: size mismatch: dim 1, 2, and 3 of matrix1 and matrix2 are unequal : ' &
1544 0 : // '(',m1,',',n1,',',j1,',',i1,') and (',m2,',',n2,',',j2,',',i2,')'
1545 0 : STOP
1546 : end if
1547 :
1548 : ! save mat1
1549 1 : call move_alloc(mat1, tmp)
1550 :
1551 6 : allocate(mat1(m1,n1,j1,i1 + i2))
1552 23 : mat1(:,:,:,1:i1) = tmp(:,:,:,1:i1)
1553 23 : mat1(:,:,:,i1+1_i4:i1+i2) = mat2(:,:,:,1:i2)
1554 : end if
1555 :
1556 : else
1557 :
1558 6 : allocate(mat1(m2,n2,j2,i2))
1559 13 : mat1 = mat2
1560 :
1561 : end if
1562 :
1563 11 : END SUBROUTINE append_lgt_4d
1564 :
1565 : ! ------------------------------------------------------------------
1566 :
1567 2 : SUBROUTINE paste_i4_m_s(mat1, sca2)
1568 :
1569 : implicit none
1570 :
1571 : integer(i4), dimension(:,:), allocatable, intent(inout) :: mat1
1572 : integer(i4), intent(in) :: sca2
1573 :
1574 : ! local variables
1575 : integer(i4) :: m1 ! dim1 of matrix
1576 : integer(i4) :: n1 ! dim2 of matrix
1577 2 : integer(i4), dimension(:,:), allocatable :: tmp
1578 :
1579 2 : if (allocated(mat1)) then
1580 1 : m1 = size(mat1,1) ! rows
1581 1 : n1 = size(mat1,2) ! columns
1582 1 : if (m1 /= 1_i4) then
1583 0 : print*, 'paste: scalar paste to matrix only works with one-line matrix'
1584 0 : STOP
1585 : end if
1586 : ! save mat1
1587 1 : call move_alloc(mat1, tmp)
1588 :
1589 3 : allocate(mat1(1_i4,n1+1_i4))
1590 2 : mat1(1,1:n1) = tmp(1,1:n1)
1591 1 : mat1(1,n1+1_i4) = sca2
1592 : else
1593 1 : allocate(mat1(1_i4,1_i4))
1594 1 : mat1(1,1) = sca2
1595 : end if
1596 :
1597 8 : END SUBROUTINE paste_i4_m_s
1598 :
1599 4 : SUBROUTINE paste_i4_m_v(mat1, vec2, fill_value)
1600 :
1601 : implicit none
1602 :
1603 : integer(i4), dimension(:,:), allocatable, intent(inout) :: mat1
1604 : integer(i4), dimension(:), intent(in) :: vec2
1605 : integer(i4), optional, intent(in) :: fill_value
1606 :
1607 : ! local variables
1608 : integer(i4) :: m1, m2 ! dim1 of matrixes
1609 : integer(i4) :: n1, n2 ! dim2 of matrixes
1610 4 : integer(i4), dimension(:,:), allocatable :: tmp
1611 :
1612 4 : m2 = size(vec2,1) ! rows
1613 4 : n2 = 1_i4 ! columns
1614 :
1615 4 : if (allocated(mat1)) then
1616 3 : m1 = size(mat1,1) ! rows
1617 3 : n1 = size(mat1,2) ! columns
1618 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
1619 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1620 0 : STOP
1621 : end if
1622 : ! save mat1
1623 3 : call move_alloc(mat1, tmp)
1624 :
1625 3 : if ( m1 == m2 ) then
1626 4 : allocate(mat1(m1,n1+n2))
1627 3 : mat1(1:m1,1:n1) = tmp(:,1:n1)
1628 2 : mat1(1:m2,n1+n2) = vec2(1:m2)
1629 : end if
1630 :
1631 3 : if ( m1 > m2 ) then
1632 4 : allocate(mat1(m1,n1+n2))
1633 13 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
1634 3 : mat1( 1:m2,n1+n2) = vec2(1:m2)
1635 2 : mat1(m2+1:m1,n1+n2) = fill_value
1636 : end if
1637 :
1638 3 : if ( m1 < m2 ) then
1639 4 : allocate(mat1(m2,n1+n2))
1640 5 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
1641 7 : mat1(m1+1:m2,1:n1) = fill_value
1642 4 : mat1( 1:m2,n1+n2) = vec2(1:m2)
1643 : end if
1644 :
1645 : else
1646 1 : n1 = 0_i4
1647 1 : m1 = m2
1648 :
1649 3 : allocate(mat1(m2,n2))
1650 3 : mat1(1:m2,n1+n2) = vec2(1:m2)
1651 : end if
1652 :
1653 6 : END SUBROUTINE paste_i4_m_v
1654 :
1655 5 : SUBROUTINE paste_i4_m_m(mat1, mat2, fill_value)
1656 :
1657 : implicit none
1658 :
1659 : integer(i4), dimension(:,:), allocatable, intent(inout) :: mat1
1660 : integer(i4), dimension(:,:), intent(in) :: mat2
1661 : integer(i4), optional, intent(in) :: fill_value
1662 :
1663 : ! local variables
1664 : integer(i4) :: m1, m2 ! dim1 of matrixes
1665 : integer(i4) :: n1, n2 ! dim2 of matrixes
1666 5 : integer(i4), dimension(:,:), allocatable :: tmp
1667 :
1668 5 : m2 = size(mat2,1) ! rows
1669 5 : n2 = size(mat2,2) ! columns
1670 :
1671 5 : if (allocated(mat1)) then
1672 4 : m1 = size(mat1,1) ! rows
1673 4 : n1 = size(mat1,2) ! columns
1674 4 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
1675 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1676 0 : STOP
1677 : end if
1678 : ! save mat1
1679 4 : call move_alloc(mat1, tmp)
1680 :
1681 4 : if ( m1 == m2 ) then
1682 8 : allocate(mat1(m1,n1+n2))
1683 23 : mat1(:,1:n1) = tmp(:,1:n1)
1684 17 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1685 : end if
1686 :
1687 4 : if ( m1 > m2 ) then
1688 4 : allocate(mat1(m1,n1+n2))
1689 7 : mat1( : ,1:n1) = tmp(:,1:n1)
1690 7 : mat1( 1:m2,n1+1_i4:n1+n2) = mat2(:,1:n2)
1691 7 : mat1(m2+1:m1,n1+1_i4:n1+n2) = fill_value
1692 : end if
1693 :
1694 4 : if ( m1 < m2 ) then
1695 4 : allocate(mat1(m2,n1+n2))
1696 16 : mat1( 1:m1, 1:n1 ) = tmp(:,1:n1)
1697 11 : mat1(m1+1:m2, 1:n1 ) = fill_value
1698 13 : mat1( : ,n1+1_i4:n1+n2) = mat2(:,1:n2)
1699 : end if
1700 :
1701 : else
1702 1 : n1 = 0_i4
1703 1 : m1 = m2
1704 :
1705 4 : allocate(mat1(m2,n2))
1706 8 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1707 : end if
1708 :
1709 9 : END SUBROUTINE paste_i4_m_m
1710 :
1711 3 : SUBROUTINE paste_i8_m_s(mat1, sca2)
1712 :
1713 : implicit none
1714 :
1715 : integer(i8), dimension(:,:), allocatable, intent(inout) :: mat1
1716 : integer(i8), intent(in) :: sca2
1717 :
1718 : ! local variables
1719 : integer(i4) :: m1 ! dim1 of matrix
1720 : integer(i4) :: n1 ! dim2 of matrix
1721 3 : integer(i8), dimension(:,:), allocatable :: tmp
1722 :
1723 3 : if (allocated(mat1)) then
1724 2 : m1 = size(mat1,1) ! rows
1725 2 : n1 = size(mat1,2) ! columns
1726 2 : if (m1 /= 1_i4) then
1727 0 : print*, 'paste: scalar paste to matrix only works with one-line matrix'
1728 0 : STOP
1729 : end if
1730 : ! save mat1
1731 2 : call move_alloc(mat1, tmp)
1732 :
1733 6 : allocate(mat1(1_i4,n1+1_i4))
1734 6 : mat1(1,1:n1) = tmp(1,1:n1)
1735 2 : mat1(1,n1+1_i4) = sca2
1736 : else
1737 1 : allocate(mat1(1_i4,1_i4))
1738 1 : mat1(1,1) = sca2
1739 : end if
1740 :
1741 8 : END SUBROUTINE paste_i8_m_s
1742 :
1743 4 : SUBROUTINE paste_i8_m_v(mat1, vec2, fill_value)
1744 :
1745 : implicit none
1746 :
1747 : integer(i8), dimension(:,:), allocatable, intent(inout) :: mat1
1748 : integer(i8), dimension(:), intent(in) :: vec2
1749 : integer(i8), optional, intent(in) :: fill_value
1750 :
1751 : ! local variables
1752 : integer(i4) :: m1, m2 ! dim1 of matrixes
1753 : integer(i4) :: n1, n2 ! dim2 of matrixes
1754 4 : integer(i8), dimension(:,:), allocatable :: tmp
1755 :
1756 4 : m2 = size(vec2,1) ! rows
1757 4 : n2 = 1_i4 ! columns
1758 :
1759 4 : if (allocated(mat1)) then
1760 3 : m1 = size(mat1,1) ! rows
1761 3 : n1 = size(mat1,2) ! columns
1762 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
1763 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1764 0 : STOP
1765 : end if
1766 : ! save mat1
1767 3 : call move_alloc(mat1, tmp)
1768 :
1769 3 : if ( m1 == m2 ) then
1770 4 : allocate(mat1(m1,n1+n2))
1771 3 : mat1(1:m1,1:n1) = tmp(:,1:n1)
1772 2 : mat1(1:m2,n1+n2) = vec2(1:m2)
1773 : end if
1774 :
1775 3 : if ( m1 > m2 ) then
1776 4 : allocate(mat1(m1,n1+n2))
1777 13 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
1778 3 : mat1( 1:m2,n1+n2) = vec2(1:m2)
1779 2 : mat1(m2+1:m1,n1+n2) = fill_value
1780 : end if
1781 :
1782 3 : if ( m1 < m2 ) then
1783 4 : allocate(mat1(m2,n1+n2))
1784 5 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
1785 7 : mat1(m1+1:m2,1:n1) = fill_value
1786 4 : mat1( 1:m2,n1+n2) = vec2(1:m2)
1787 : end if
1788 :
1789 : else
1790 1 : n1 = 0_i4
1791 1 : m1 = m2
1792 :
1793 3 : allocate(mat1(m2,n2))
1794 3 : mat1(1:m2,n1+n2) = vec2(1:m2)
1795 : end if
1796 :
1797 7 : END SUBROUTINE paste_i8_m_v
1798 :
1799 4 : SUBROUTINE paste_i8_m_m(mat1, mat2, fill_value)
1800 :
1801 : implicit none
1802 :
1803 : integer(i8), dimension(:,:), allocatable, intent(inout) :: mat1
1804 : integer(i8), dimension(:,:), intent(in) :: mat2
1805 : integer(i8), optional, intent(in) :: fill_value
1806 :
1807 : ! local variables
1808 : integer(i4) :: m1, m2 ! dim1 of matrixes
1809 : integer(i4) :: n1, n2 ! dim2 of matrixes
1810 4 : integer(i8), dimension(:,:), allocatable :: tmp
1811 :
1812 4 : m2 = size(mat2,1) ! rows
1813 4 : n2 = size(mat2,2) ! columns
1814 :
1815 4 : if (allocated(mat1)) then
1816 3 : m1 = size(mat1,1) ! rows
1817 3 : n1 = size(mat1,2) ! columns
1818 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
1819 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1820 0 : STOP
1821 : end if
1822 : ! save mat1
1823 3 : call move_alloc(mat1, tmp)
1824 :
1825 3 : if ( m1 == m2 ) then
1826 4 : allocate(mat1(m1,n1+n2))
1827 7 : mat1(:,1:n1) = tmp(:,1:n1)
1828 10 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1829 : end if
1830 :
1831 3 : if ( m1 > m2 ) then
1832 4 : allocate(mat1(m1,n1+n2))
1833 7 : mat1( : ,1:n1) = tmp(:,1:n1)
1834 7 : mat1( 1:m2,n1+1_i4:n1+n2) = mat2(:,1:n2)
1835 7 : mat1(m2+1:m1,n1+1_i4:n1+n2) = fill_value
1836 : end if
1837 :
1838 3 : if ( m1 < m2 ) then
1839 4 : allocate(mat1(m2,n1+n2))
1840 16 : mat1( 1:m1, 1:n1 ) = tmp(:,1:n1)
1841 11 : mat1(m1+1:m2, 1:n1 ) = fill_value
1842 13 : mat1( : ,n1+1_i4:n1+n2) = mat2(:,1:n2)
1843 : end if
1844 :
1845 : else
1846 1 : n1 = 0_i4
1847 1 : m1 = m2
1848 :
1849 4 : allocate(mat1(m2,n2))
1850 8 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1851 : end if
1852 :
1853 8 : END SUBROUTINE paste_i8_m_m
1854 :
1855 2 : SUBROUTINE paste_sp_m_s(mat1, sca2)
1856 :
1857 : implicit none
1858 :
1859 : real(sp), dimension(:,:), allocatable, intent(inout) :: mat1
1860 : real(sp), intent(in) :: sca2
1861 :
1862 : ! local variables
1863 : integer(i4) :: m1 ! dim1 of matrix
1864 : integer(i4) :: n1 ! dim2 of matrix
1865 2 : real(sp), dimension(:,:), allocatable :: tmp
1866 :
1867 2 : if (allocated(mat1)) then
1868 1 : m1 = size(mat1,1) ! rows
1869 1 : n1 = size(mat1,2) ! columns
1870 1 : if (m1 /= 1_i4) then
1871 0 : print*, 'paste: scalar paste to matrix only works with one-line matrix'
1872 0 : STOP
1873 : end if
1874 : ! save mat1
1875 1 : call move_alloc(mat1, tmp)
1876 :
1877 3 : allocate(mat1(1_i4,n1+1_i4))
1878 2 : mat1(1,1:n1) = tmp(1,1:n1)
1879 1 : mat1(1,n1+1_i4) = sca2
1880 : else
1881 1 : allocate(mat1(1_i4,1_i4))
1882 1 : mat1(1,1) = sca2
1883 : end if
1884 :
1885 6 : END SUBROUTINE paste_sp_m_s
1886 :
1887 4 : SUBROUTINE paste_sp_m_v(mat1, vec2, fill_value)
1888 :
1889 : implicit none
1890 :
1891 : real(sp), dimension(:,:), allocatable, intent(inout) :: mat1
1892 : real(sp), dimension(:), intent(in) :: vec2
1893 : real(sp), optional, intent(in) :: fill_value
1894 :
1895 : ! local variables
1896 : integer(i4) :: m1, m2 ! dim1 of matrixes
1897 : integer(i4) :: n1, n2 ! dim2 of matrixes
1898 4 : real(sp), dimension(:,:), allocatable :: tmp
1899 :
1900 4 : m2 = size(vec2,1) ! rows
1901 4 : n2 = 1_i4 ! columns
1902 :
1903 4 : if (allocated(mat1)) then
1904 3 : m1 = size(mat1,1) ! rows
1905 3 : n1 = size(mat1,2) ! columns
1906 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
1907 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1908 0 : STOP
1909 : end if
1910 : ! save mat1
1911 3 : call move_alloc(mat1, tmp)
1912 :
1913 3 : if ( m1 == m2 ) then
1914 4 : allocate(mat1(m1,n1+n2))
1915 3 : mat1(1:m1,1:n1) = tmp(:,1:n1)
1916 2 : mat1(1:m2,n1+n2) = vec2(1:m2)
1917 : end if
1918 :
1919 3 : if ( m1 > m2 ) then
1920 4 : allocate(mat1(m1,n1+n2))
1921 13 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
1922 3 : mat1( 1:m2,n1+n2) = vec2(1:m2)
1923 2 : mat1(m2+1:m1,n1+n2) = fill_value
1924 : end if
1925 :
1926 3 : if ( m1 < m2 ) then
1927 4 : allocate(mat1(m2,n1+n2))
1928 5 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
1929 7 : mat1(m1+1:m2,1:n1) = fill_value
1930 4 : mat1( 1:m2,n1+n2) = vec2(1:m2)
1931 : end if
1932 :
1933 : else
1934 1 : n1 = 0_i4
1935 1 : m1 = m2
1936 :
1937 3 : allocate(mat1(m2,n2))
1938 3 : mat1(1:m2,n1+n2) = vec2(1:m2)
1939 : end if
1940 :
1941 6 : END SUBROUTINE paste_sp_m_v
1942 :
1943 4 : SUBROUTINE paste_sp_m_m(mat1, mat2, fill_value)
1944 :
1945 : implicit none
1946 :
1947 : real(sp), dimension(:,:), allocatable, intent(inout) :: mat1
1948 : real(sp), dimension(:,:), intent(in) :: mat2
1949 : real(sp), optional, intent(in) :: fill_value
1950 :
1951 : ! local variables
1952 : integer(i4) :: m1, m2 ! dim1 of matrixes
1953 : integer(i4) :: n1, n2 ! dim2 of matrixes
1954 4 : real(sp), dimension(:,:), allocatable :: tmp
1955 :
1956 4 : m2 = size(mat2,1) ! rows
1957 4 : n2 = size(mat2,2) ! columns
1958 :
1959 4 : if (allocated(mat1)) then
1960 3 : m1 = size(mat1,1) ! rows
1961 3 : n1 = size(mat1,2) ! columns
1962 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
1963 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
1964 0 : STOP
1965 : end if
1966 : ! save mat1
1967 3 : call move_alloc(mat1, tmp)
1968 :
1969 3 : if ( m1 == m2 ) then
1970 4 : allocate(mat1(m1,n1+n2))
1971 7 : mat1(:,1:n1) = tmp(:,1:n1)
1972 10 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1973 : end if
1974 :
1975 3 : if ( m1 > m2 ) then
1976 4 : allocate(mat1(m1,n1+n2))
1977 7 : mat1( : ,1:n1) = tmp(:,1:n1)
1978 7 : mat1( 1:m2,n1+1_i4:n1+n2) = mat2(:,1:n2)
1979 7 : mat1(m2+1:m1,n1+1_i4:n1+n2) = fill_value
1980 : end if
1981 :
1982 3 : if ( m1 < m2 ) then
1983 4 : allocate(mat1(m2,n1+n2))
1984 16 : mat1( 1:m1, 1:n1 ) = tmp(:,1:n1)
1985 11 : mat1(m1+1:m2, 1:n1 ) = fill_value
1986 13 : mat1( : ,n1+1_i4:n1+n2) = mat2(:,1:n2)
1987 : end if
1988 :
1989 :
1990 : else
1991 1 : n1 = 0_i4
1992 1 : m1 = m2
1993 :
1994 4 : allocate(mat1(m2,n2))
1995 8 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
1996 : end if
1997 :
1998 8 : END SUBROUTINE paste_sp_m_m
1999 :
2000 2 : SUBROUTINE paste_dp_m_s(mat1, sca2)
2001 :
2002 : implicit none
2003 :
2004 : real(dp), dimension(:,:), allocatable, intent(inout) :: mat1
2005 : real(dp), intent(in) :: sca2
2006 :
2007 : ! local variables
2008 : integer(i4) :: m1 ! dim1 of matrix
2009 : integer(i4) :: n1 ! dim2 of matrix
2010 2 : real(dp), dimension(:,:), allocatable :: tmp
2011 :
2012 2 : if (allocated(mat1)) then
2013 1 : m1 = size(mat1,1) ! rows
2014 1 : n1 = size(mat1,2) ! columns
2015 1 : if (m1 /= 1_i4) then
2016 0 : print*, 'paste: scalar paste to matrix only works with one-line matrix'
2017 0 : STOP
2018 : end if
2019 : ! save mat1
2020 1 : call move_alloc(mat1, tmp)
2021 :
2022 3 : allocate(mat1(1_i4,n1+1_i4))
2023 2 : mat1(1,1:n1) = tmp(1,1:n1)
2024 1 : mat1(1,n1+1_i4) = sca2
2025 : else
2026 1 : allocate(mat1(1_i4,1_i4))
2027 1 : mat1(1,1) = sca2
2028 : end if
2029 :
2030 6 : END SUBROUTINE paste_dp_m_s
2031 :
2032 4 : SUBROUTINE paste_dp_m_v(mat1, vec2, fill_value)
2033 :
2034 : implicit none
2035 :
2036 : real(dp), dimension(:,:), allocatable, intent(inout) :: mat1
2037 : real(dp), dimension(:), intent(in) :: vec2
2038 : real(dp), optional, intent(in) :: fill_value
2039 :
2040 : ! local variables
2041 : integer(i4) :: m1, m2 ! dim1 of matrixes
2042 : integer(i4) :: n1, n2 ! dim2 of matrixes
2043 4 : real(dp), dimension(:,:), allocatable :: tmp
2044 :
2045 4 : m2 = size(vec2,1) ! rows
2046 4 : n2 = 1_i4 ! columns
2047 :
2048 4 : if (allocated(mat1)) then
2049 3 : m1 = size(mat1,1) ! rows
2050 3 : n1 = size(mat1,2) ! columns
2051 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
2052 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
2053 0 : STOP
2054 : end if
2055 : ! save mat1
2056 3 : call move_alloc(mat1, tmp)
2057 :
2058 3 : if ( m1 == m2 ) then
2059 4 : allocate(mat1(m1,n1+n2))
2060 3 : mat1(1:m1,1:n1) = tmp(:,1:n1)
2061 2 : mat1(1:m2,n1+n2) = vec2(1:m2)
2062 : end if
2063 :
2064 3 : if ( m1 > m2 ) then
2065 4 : allocate(mat1(m1,n1+n2))
2066 13 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
2067 3 : mat1( 1:m2,n1+n2) = vec2(1:m2)
2068 2 : mat1(m2+1:m1,n1+n2) = fill_value
2069 : end if
2070 :
2071 3 : if ( m1 < m2 ) then
2072 4 : allocate(mat1(m2,n1+n2))
2073 5 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
2074 7 : mat1(m1+1:m2,1:n1) = fill_value
2075 4 : mat1( 1:m2,n1+n2) = vec2(1:m2)
2076 : end if
2077 :
2078 : else
2079 1 : n1 = 0_i4
2080 1 : m1 = m2
2081 :
2082 3 : allocate(mat1(m2,n2))
2083 3 : mat1(1:m2,n1+n2) = vec2(1:m2)
2084 : end if
2085 :
2086 6 : END SUBROUTINE paste_dp_m_v
2087 :
2088 5 : SUBROUTINE paste_dp_m_m(mat1, mat2, fill_value)
2089 :
2090 : implicit none
2091 :
2092 : real(dp), dimension(:,:), allocatable, intent(inout) :: mat1
2093 : real(dp), dimension(:,:), intent(in) :: mat2
2094 : real(dp), optional, intent(in) :: fill_value
2095 :
2096 : ! local variables
2097 : integer(i4) :: m1, m2 ! dim1 of matrixes
2098 : integer(i4) :: n1, n2 ! dim2 of matrixes
2099 5 : real(dp), dimension(:,:), allocatable :: tmp
2100 :
2101 5 : m2 = size(mat2,1) ! rows
2102 5 : n2 = size(mat2,2) ! columns
2103 :
2104 5 : if (allocated(mat1)) then
2105 4 : m1 = size(mat1,1) ! rows
2106 4 : n1 = size(mat1,2) ! columns
2107 4 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
2108 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
2109 0 : STOP
2110 : end if
2111 : ! save mat1
2112 4 : call move_alloc(mat1, tmp)
2113 :
2114 4 : if ( m1 == m2 ) then
2115 8 : allocate(mat1(m1,n1+n2))
2116 23 : mat1(:,1:n1) = tmp(:,1:n1)
2117 17 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
2118 : end if
2119 :
2120 4 : if ( m1 > m2 ) then
2121 4 : allocate(mat1(m1,n1+n2))
2122 7 : mat1( : ,1:n1) = tmp(:,1:n1)
2123 7 : mat1( 1:m2,n1+1_i4:n1+n2) = mat2(:,1:n2)
2124 7 : mat1(m2+1:m1,n1+1_i4:n1+n2) = fill_value
2125 : end if
2126 :
2127 4 : if ( m1 < m2 ) then
2128 4 : allocate(mat1(m2,n1+n2))
2129 16 : mat1( 1:m1, 1:n1 ) = tmp(:,1:n1)
2130 11 : mat1(m1+1:m2, 1:n1 ) = fill_value
2131 13 : mat1( : ,n1+1_i4:n1+n2) = mat2(:,1:n2)
2132 : end if
2133 :
2134 : else
2135 1 : n1 = 0_i4
2136 1 : m1 = m2
2137 :
2138 4 : allocate(mat1(m2,n2))
2139 8 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
2140 : end if
2141 :
2142 9 : END SUBROUTINE paste_dp_m_m
2143 :
2144 3 : SUBROUTINE paste_dp_3d(mat1, mat2)
2145 :
2146 : real(dp), dimension(:,:,:), allocatable, intent(inout) :: mat1
2147 : real(dp), dimension(:,:,:), intent(in) :: mat2
2148 :
2149 : ! local variables
2150 : integer(i4) :: m1, m2 ! dim1 of matrixes
2151 : integer(i4) :: n1, n2 ! dim2 of matrixes
2152 : integer(i4) :: o1, o2 ! dim2 of matrixes
2153 3 : real(dp), dimension(:,:,:), allocatable :: tmp
2154 :
2155 3 : m2 = size(mat2,1) ! rows
2156 3 : n2 = size(mat2,2) ! columns
2157 3 : o2 = size(mat2,3) ! 3rd
2158 :
2159 3 : if (allocated(mat1)) then
2160 2 : m1 = size(mat1,1) ! rows
2161 2 : n1 = size(mat1,2) ! columns
2162 2 : o1 = size(mat1,3) ! columns
2163 2 : if ( (m1 /= m2)) then
2164 0 : print*, 'paste: rows of array1 and array2 are unequal : (',m1,') and (',m2,')'
2165 0 : STOP
2166 2 : else if ( (n1 /= n2)) then
2167 0 : print*, 'paste: columns of array1 and array2 are unequal : (',n1,') and (',n2,')'
2168 0 : STOP
2169 : end if
2170 : ! save mat1
2171 2 : call move_alloc(mat1, tmp)
2172 :
2173 10 : allocate(mat1(m1,n1,o1+o2))
2174 37 : mat1(:,:,1:o1) = tmp(:,:,:)
2175 23 : mat1(:,:,o1+1_i4:o1+o2) = mat2(:,:,:)
2176 :
2177 : else
2178 5 : allocate(mat1(m2,n2,o2))
2179 16 : mat1(:,:,:) = mat2(:,:,:)
2180 : end if
2181 :
2182 8 : END SUBROUTINE paste_dp_3d
2183 :
2184 3 : SUBROUTINE paste_dp_4d(mat1, mat2)
2185 :
2186 : real(dp), dimension(:,:,:,:), allocatable, intent(inout) :: mat1
2187 : real(dp), dimension(:,:,:,:), intent(in) :: mat2
2188 :
2189 : ! local variables
2190 : integer(i4) :: m1, m2 ! dim1 of matrixes
2191 : integer(i4) :: n1, n2 ! dim2 of matrixes
2192 : integer(i4) :: o1, o2 ! dim3 of matrixes
2193 : integer(i4) :: p1, p2 ! dim4 of matrixes
2194 3 : real(dp), dimension(:,:,:,:), allocatable :: tmp
2195 :
2196 3 : m2 = size(mat2,1) ! rows
2197 3 : n2 = size(mat2,2) ! columns
2198 3 : o2 = size(mat2,3) ! 3rd
2199 3 : p2 = size(mat2,4) ! 4th
2200 :
2201 3 : if (allocated(mat1)) then
2202 2 : m1 = size(mat1,1) ! rows
2203 2 : n1 = size(mat1,2) ! columns
2204 2 : o1 = size(mat1,3) ! 3rd
2205 2 : p1 = size(mat1,4) ! 4th
2206 2 : if ( (m1 /= m2)) then
2207 0 : print*, 'paste: rows of array1 and array2 are unequal : (',m1,') and (',m2,')'
2208 0 : STOP
2209 2 : else if ( (n1 /= n2)) then
2210 0 : print*, 'paste: columns of array1 and array2 are unequal : (',n1,') and (',n2,')'
2211 0 : STOP
2212 2 : else if ( (o1 /= o2)) then
2213 0 : print*, 'paste: columns of array1 and array2 are unequal : (',o1,') and (',o2,')'
2214 0 : STOP
2215 : end if
2216 : ! save mat1
2217 2 : call move_alloc(mat1, tmp)
2218 :
2219 12 : allocate(mat1(m1,n1,o1,p1+p2))
2220 77 : mat1(:,:,:,1:p1) = tmp(:,:,:,:)
2221 47 : mat1(:,:,:,p1+1_i4:p1+p2) = mat2(:,:,:,:)
2222 :
2223 : else
2224 6 : allocate(mat1(m2,n2,o2,p2))
2225 32 : mat1(:,:,:,:) = mat2(:,:,:,:)
2226 : end if
2227 :
2228 6 : END SUBROUTINE paste_dp_4d
2229 :
2230 3 : SUBROUTINE paste_i4_3d(mat1, mat2)
2231 :
2232 : integer(i4), dimension(:,:,:), allocatable, intent(inout) :: mat1
2233 : integer(i4), dimension(:,:,:), intent(in) :: mat2
2234 :
2235 : ! local variables
2236 : integer(i4) :: m1, m2 ! dim1 of matrixes
2237 : integer(i4) :: n1, n2 ! dim2 of matrixes
2238 : integer(i4) :: o1, o2 ! dim2 of matrixes
2239 3 : integer(i4), dimension(:,:,:), allocatable :: tmp
2240 :
2241 3 : m2 = size(mat2,1) ! rows
2242 3 : n2 = size(mat2,2) ! columns
2243 3 : o2 = size(mat2,3) ! 3rd
2244 :
2245 3 : if (allocated(mat1)) then
2246 2 : m1 = size(mat1,1) ! rows
2247 2 : n1 = size(mat1,2) ! columns
2248 2 : o1 = size(mat1,3) ! columns
2249 2 : if ( (m1 /= m2)) then
2250 0 : print*, 'paste: rows of array1 and array2 are unequal : (',m1,') and (',m2,')'
2251 0 : STOP
2252 2 : else if ( (n1 /= n2)) then
2253 0 : print*, 'paste: columns of array1 and array2 are unequal : (',n1,') and (',n2,')'
2254 0 : STOP
2255 : end if
2256 : ! save mat1
2257 2 : call move_alloc(mat1, tmp)
2258 :
2259 10 : allocate(mat1(m1,n1,o1+o2))
2260 37 : mat1(:,:,1:o1) = tmp(:,:,:)
2261 23 : mat1(:,:,o1+1_i4:o1+o2) = mat2(:,:,:)
2262 :
2263 : else
2264 5 : allocate(mat1(m2,n2,o2))
2265 16 : mat1(:,:,:) = mat2(:,:,:)
2266 : end if
2267 :
2268 6 : END SUBROUTINE paste_i4_3d
2269 :
2270 3 : SUBROUTINE paste_i4_4d(mat1, mat2)
2271 :
2272 : integer(i4), dimension(:,:,:,:), allocatable, intent(inout) :: mat1
2273 : integer(i4), dimension(:,:,:,:), intent(in) :: mat2
2274 :
2275 : ! local variables
2276 : integer(i4) :: m1, m2 ! dim1 of matrixes
2277 : integer(i4) :: n1, n2 ! dim2 of matrixes
2278 : integer(i4) :: o1, o2 ! dim3 of matrixes
2279 : integer(i4) :: p1, p2 ! dim4 of matrixes
2280 3 : integer(i4), dimension(:,:,:,:), allocatable :: tmp
2281 :
2282 3 : m2 = size(mat2,1) ! rows
2283 3 : n2 = size(mat2,2) ! columns
2284 3 : o2 = size(mat2,3) ! 3rd
2285 3 : p2 = size(mat2,4) ! 4th
2286 :
2287 3 : if (allocated(mat1)) then
2288 2 : m1 = size(mat1,1) ! rows
2289 2 : n1 = size(mat1,2) ! columns
2290 2 : o1 = size(mat1,3) ! 3rd
2291 2 : p1 = size(mat1,4) ! 4th
2292 2 : if ( (m1 /= m2)) then
2293 0 : print*, 'paste: rows of array1 and array2 are unequal : (',m1,') and (',m2,')'
2294 0 : STOP
2295 2 : else if ( (n1 /= n2)) then
2296 0 : print*, 'paste: columns of array1 and array2 are unequal : (',n1,') and (',n2,')'
2297 0 : STOP
2298 2 : else if ( (o1 /= o2)) then
2299 0 : print*, 'paste: columns of array1 and array2 are unequal : (',o1,') and (',o2,')'
2300 0 : STOP
2301 : end if
2302 : ! save mat1
2303 2 : call move_alloc(mat1, tmp)
2304 :
2305 12 : allocate(mat1(m1,n1,o1,p1+p2))
2306 77 : mat1(:,:,:,1:p1) = tmp(:,:,:,:)
2307 47 : mat1(:,:,:,p1+1_i4:p1+p2) = mat2(:,:,:,:)
2308 :
2309 : else
2310 6 : allocate(mat1(m2,n2,o2,p2))
2311 32 : mat1(:,:,:,:) = mat2(:,:,:,:)
2312 : end if
2313 :
2314 6 : END SUBROUTINE paste_i4_4d
2315 :
2316 2 : SUBROUTINE paste_char_m_s(mat1, sca2)
2317 :
2318 : implicit none
2319 :
2320 : character(len=*), dimension(:,:), allocatable, intent(inout) :: mat1
2321 : character(len=*), intent(in) :: sca2
2322 :
2323 : ! local variables
2324 : integer(i4) :: m1 ! dim1 of matrix
2325 : integer(i4) :: n1 ! dim2 of matrix
2326 2 : character(len(mat1)), dimension(:,:), allocatable :: tmp
2327 :
2328 2 : if (allocated(mat1)) then
2329 1 : m1 = size(mat1,1) ! rows
2330 1 : n1 = size(mat1,2) ! columns
2331 1 : if (m1 /= 1_i4) then
2332 0 : print*, 'paste: scalar paste to matrix only works with one-line matrix'
2333 0 : STOP
2334 : end if
2335 : ! save mat1
2336 1 : call move_alloc(mat1, tmp)
2337 :
2338 3 : allocate(mat1(1_i4,n1+1_i4))
2339 2 : mat1(1,1:n1) = tmp(1,1:n1)
2340 1 : mat1(1,n1+1_i4) = sca2
2341 : else
2342 1 : allocate(mat1(1_i4,1_i4))
2343 1 : mat1(1,1) = sca2
2344 : end if
2345 :
2346 5 : END SUBROUTINE paste_char_m_s
2347 :
2348 4 : SUBROUTINE paste_char_m_v(mat1, vec2, fill_value)
2349 :
2350 : implicit none
2351 :
2352 : character(len=*), dimension(:,:), allocatable, intent(inout) :: mat1
2353 : character(len=*), dimension(:), intent(in) :: vec2
2354 : character(len=*), optional, intent(in) :: fill_value
2355 :
2356 : ! local variables
2357 : integer(i4) :: m1, m2 ! dim1 of matrixes
2358 : integer(i4) :: n1, n2 ! dim2 of matrixes
2359 4 : character(len(mat1)), dimension(:,:), allocatable :: tmp
2360 :
2361 4 : m2 = size(vec2,1) ! rows
2362 4 : n2 = 1_i4 ! columns
2363 :
2364 4 : if (allocated(mat1)) then
2365 3 : m1 = size(mat1,1) ! rows
2366 3 : n1 = size(mat1,2) ! columns
2367 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
2368 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
2369 0 : STOP
2370 : end if
2371 : ! save mat1
2372 3 : call move_alloc(mat1, tmp)
2373 :
2374 3 : if ( m1 == m2 ) then
2375 4 : allocate(mat1(m1,n1+n2))
2376 3 : mat1(1:m1,1:n1) = tmp(:,1:n1)
2377 2 : mat1(1:m2,n1+n2) = vec2(1:m2)
2378 : end if
2379 :
2380 3 : if ( m1 > m2 ) then
2381 4 : allocate(mat1(m1,n1+n2))
2382 13 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
2383 3 : mat1( 1:m2,n1+n2) = vec2(1:m2)
2384 2 : mat1(m2+1:m1,n1+n2) = fill_value
2385 : end if
2386 :
2387 3 : if ( m1 < m2 ) then
2388 4 : allocate(mat1(m2,n1+n2))
2389 5 : mat1( 1:m1,1:n1) = tmp(:,1:n1)
2390 8 : mat1(m1+1:m2,1:n1) = fill_value
2391 4 : mat1( 1:m2,n1+n2) = vec2(1:m2)
2392 : end if
2393 :
2394 : else
2395 1 : n1 = 0_i4
2396 1 : m1 = m2
2397 :
2398 3 : allocate(mat1(m2,n2))
2399 3 : mat1(1:m2,n1+n2) = vec2(1:m2)
2400 : end if
2401 :
2402 6 : END SUBROUTINE paste_char_m_v
2403 :
2404 4 : SUBROUTINE paste_char_m_m(mat1, mat2, fill_value)
2405 :
2406 : implicit none
2407 :
2408 : character(len=*), dimension(:,:), allocatable, intent(inout) :: mat1
2409 : character(len=*), dimension(:,:), intent(in) :: mat2
2410 : character(len=*), optional, intent(in) :: fill_value
2411 :
2412 : ! local variables
2413 : integer(i4) :: m1, m2 ! dim1 of matrixes
2414 : integer(i4) :: n1, n2 ! dim2 of matrixes
2415 4 : character(len(mat1)), dimension(:,:), allocatable :: tmp
2416 :
2417 4 : m2 = size(mat2,1) ! rows
2418 4 : n2 = size(mat2,2) ! columns
2419 :
2420 4 : if (allocated(mat1)) then
2421 3 : m1 = size(mat1,1) ! rows
2422 3 : n1 = size(mat1,2) ! columns
2423 3 : if ( (m1 /= m2) .and. .not. present( fill_value ) ) then
2424 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
2425 0 : STOP
2426 : end if
2427 : ! save mat1
2428 3 : call move_alloc(mat1, tmp)
2429 :
2430 3 : if ( m1 == m2 ) then
2431 4 : allocate(mat1(m1,n1+n2))
2432 7 : mat1(:,1:n1) = tmp(:,1:n1)
2433 10 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
2434 : end if
2435 :
2436 3 : if ( m1 > m2 ) then
2437 4 : allocate(mat1(m1,n1+n2))
2438 7 : mat1( : ,1:n1) = tmp(:,1:n1)
2439 7 : mat1( 1:m2,n1+1_i4:n1+n2) = mat2(:,1:n2)
2440 8 : mat1(m2+1:m1,n1+1_i4:n1+n2) = fill_value
2441 : end if
2442 :
2443 3 : if ( m1 < m2 ) then
2444 4 : allocate(mat1(m2,n1+n2))
2445 16 : mat1( 1:m1, 1:n1 ) = tmp(:,1:n1)
2446 12 : mat1(m1+1:m2, 1:n1 ) = fill_value
2447 13 : mat1( : ,n1+1_i4:n1+n2) = mat2(:,1:n2)
2448 : end if
2449 :
2450 : else
2451 1 : n1 = 0_i4
2452 1 : m1 = m2
2453 :
2454 4 : allocate(mat1(m2,n2))
2455 8 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
2456 : end if
2457 :
2458 8 : END SUBROUTINE paste_char_m_m
2459 :
2460 2 : SUBROUTINE paste_lgt_m_s(mat1, sca2)
2461 :
2462 : implicit none
2463 :
2464 : logical, dimension(:,:), allocatable, intent(inout) :: mat1
2465 : logical, intent(in) :: sca2
2466 :
2467 : ! local variables
2468 : integer(i4) :: m1 ! dim1 of matrix
2469 : integer(i4) :: n1 ! dim2 of matrix
2470 2 : logical, dimension(:,:), allocatable :: tmp
2471 :
2472 2 : if (allocated(mat1)) then
2473 1 : m1 = size(mat1,1) ! rows
2474 1 : n1 = size(mat1,2) ! columns
2475 1 : if (m1 /= 1_i4) then
2476 0 : print*, 'paste: scalar paste to matrix only works with one-line matrix'
2477 0 : STOP
2478 : end if
2479 : ! save mat1
2480 1 : call move_alloc(mat1, tmp)
2481 :
2482 3 : allocate(mat1(1_i4,n1+1_i4))
2483 2 : mat1(1,1:n1) = tmp(1,1:n1)
2484 1 : mat1(1,n1+1_i4) = sca2
2485 : else
2486 1 : allocate(mat1(1_i4,1_i4))
2487 1 : mat1(1,1) = sca2
2488 : end if
2489 :
2490 6 : END SUBROUTINE paste_lgt_m_s
2491 :
2492 2 : SUBROUTINE paste_lgt_m_v(mat1, vec2)
2493 :
2494 : implicit none
2495 :
2496 : logical, dimension(:,:), allocatable, intent(inout) :: mat1
2497 : logical, dimension(:), intent(in) :: vec2
2498 :
2499 : ! local variables
2500 : integer(i4) :: m1, m2 ! dim1 of matrixes
2501 : integer(i4) :: n1, n2 ! dim2 of matrixes
2502 2 : logical, dimension(:,:), allocatable :: tmp
2503 :
2504 2 : m2 = size(vec2,1) ! rows
2505 2 : n2 = 1_i4 ! columns
2506 :
2507 2 : if (allocated(mat1)) then
2508 1 : m1 = size(mat1,1) ! rows
2509 1 : n1 = size(mat1,2) ! columns
2510 1 : if (m1 /= m2) then
2511 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
2512 0 : STOP
2513 : end if
2514 : ! save mat1
2515 1 : call move_alloc(mat1, tmp)
2516 :
2517 4 : allocate(mat1(m1,n1+n2))
2518 3 : mat1(:,1:n1) = tmp(:,1:n1)
2519 2 : mat1(1:m2,n1+n2) = vec2(1:m2)
2520 : else
2521 1 : n1 = 0_i4
2522 1 : m1 = m2
2523 :
2524 3 : allocate(mat1(m2,n2))
2525 3 : mat1(1:m2,n1+n2) = vec2(1:m2)
2526 : end if
2527 :
2528 4 : END SUBROUTINE paste_lgt_m_v
2529 :
2530 2 : SUBROUTINE paste_lgt_m_m(mat1, mat2)
2531 :
2532 : implicit none
2533 :
2534 : logical, dimension(:,:), allocatable, intent(inout) :: mat1
2535 : logical, dimension(:,:), intent(in) :: mat2
2536 :
2537 : ! local variables
2538 : integer(i4) :: m1, m2 ! dim1 of matrixes
2539 : integer(i4) :: n1, n2 ! dim2 of matrixes
2540 2 : logical, dimension(:,:), allocatable :: tmp
2541 :
2542 2 : m2 = size(mat2,1) ! rows
2543 2 : n2 = size(mat2,2) ! columns
2544 :
2545 2 : if (allocated(mat1)) then
2546 1 : m1 = size(mat1,1) ! rows
2547 1 : n1 = size(mat1,2) ! columns
2548 1 : if (m1 /= m2) then
2549 0 : print*, 'paste: rows of matrix1 and matrix2 are unequal : (',m1,',',n1,') and (',m2,',',n2,')'
2550 0 : STOP
2551 : end if
2552 : ! save mat1
2553 1 : call move_alloc(mat1, tmp)
2554 :
2555 4 : allocate(mat1(m1,n1+n2))
2556 7 : mat1(:,1:n1) = tmp(:,1:n1)
2557 10 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
2558 : else
2559 1 : n1 = 0_i4
2560 1 : m1 = m2
2561 :
2562 4 : allocate(mat1(m2,n2))
2563 8 : mat1(:,n1+1_i4:n1+n2) = mat2(:,1:n2)
2564 : end if
2565 :
2566 4 : END SUBROUTINE paste_lgt_m_m
2567 :
2568 1 : subroutine add_nodata_slice_dp_2d(array, nAdd, noDataValue)
2569 : real(dp), dimension(:, :), intent(inout), allocatable :: array
2570 : integer(i4), intent(in) :: nAdd
2571 : real(dp), intent(in) :: noDataValue
2572 :
2573 7 : real(dp), dimension(size(array, 1), nAdd) :: dummy
2574 :
2575 1 : if (nAdd > 0_i4) then
2576 7 : dummy = noDataValue
2577 1 : call paste(array, dummy)
2578 : end if
2579 :
2580 2 : end subroutine add_nodata_slice_dp_2d
2581 :
2582 1 : subroutine add_nodata_slice_dp_3d(array, nAdd, noDataValue)
2583 : real(dp), dimension(:, :, :), intent(inout), allocatable :: array
2584 : integer(i4), intent(in) :: nAdd
2585 : real(dp), intent(in) :: noDataValue
2586 :
2587 15 : real(dp), dimension(size(array, 1), size(array, 2), nAdd) :: dummy
2588 :
2589 1 : if (nAdd > 0_i4) then
2590 15 : dummy = noDataValue
2591 1 : call paste(array, dummy)
2592 : end if
2593 :
2594 1 : end subroutine add_nodata_slice_dp_3d
2595 :
2596 1 : subroutine add_nodata_slice_dp_4d(array, nAdd, noDataValue)
2597 : real(dp), dimension(:, :, :, :), intent(inout), allocatable :: array
2598 : integer(i4), intent(in) :: nAdd
2599 : real(dp), intent(in) :: noDataValue
2600 :
2601 31 : real(dp), dimension(size(array, 1), size(array, 2), size(array, 3), nAdd):: dummy
2602 :
2603 1 : if (nAdd > 0_i4) then
2604 31 : dummy = noDataValue
2605 1 : call paste(array, dummy)
2606 : end if
2607 :
2608 1 : end subroutine add_nodata_slice_dp_4d
2609 :
2610 1 : subroutine add_nodata_slice_i4_2d(array, nAdd, noDataValue)
2611 : integer(i4), dimension(:, :), intent(inout), allocatable :: array
2612 : integer(i4), intent(in) :: nAdd
2613 : integer(i4), intent(in) :: noDataValue
2614 :
2615 1 : integer(i4), dimension(size(array, 1), nAdd) :: dummy
2616 :
2617 1 : if (nAdd > 0_i4) then
2618 7 : dummy = noDataValue
2619 1 : call paste(array, dummy)
2620 : end if
2621 :
2622 1 : end subroutine add_nodata_slice_i4_2d
2623 :
2624 1 : subroutine add_nodata_slice_i4_3d(array, nAdd, noDataValue)
2625 : integer(i4), dimension(:, :, :), intent(inout), allocatable :: array
2626 : integer(i4), intent(in) :: nAdd
2627 : integer(i4), intent(in) :: noDataValue
2628 :
2629 1 : integer(i4), dimension(size(array, 1), size(array, 2), nAdd) :: dummy
2630 :
2631 1 : if (nAdd > 0_i4) then
2632 15 : dummy = noDataValue
2633 1 : call paste(array, dummy)
2634 : end if
2635 :
2636 1 : end subroutine add_nodata_slice_i4_3d
2637 :
2638 1 : subroutine add_nodata_slice_i4_4d(array, nAdd, noDataValue)
2639 : integer(i4), dimension(:, :, :, :), intent(inout), allocatable :: array
2640 : integer(i4), intent(in) :: nAdd
2641 : integer(i4), intent(in) :: noDataValue
2642 :
2643 1 : integer(i4), dimension(size(array, 1), size(array, 2), size(array, 3), nAdd):: dummy
2644 :
2645 1 : if (nAdd > 0_i4) then
2646 31 : dummy = noDataValue
2647 1 : call paste(array, dummy)
2648 : end if
2649 :
2650 1 : end subroutine add_nodata_slice_i4_4d
2651 :
2652 : END MODULE mo_append
|