1module bMatrix_mod
2 ! MODULE bMatrix_mod (prefix='bmat' category='2. B and R matrices')
3 !
4 !:Purpose: A higher-level module that takes care of calling subroutines
5 ! in the lower-level modules for each specific type of B matrix:
6 ! * bMatrixHI/lamBmatrixHI
7 ! * bMatrixEnsemble
8 ! * bMatrixDiff
9 ! * bMatrixChem
10 !
11 !:Comments:
12 ! - Considerations for ensemble-based and regional static
13 ! covariances for constituents are not yet included.
14 !
15 use midasMpi_mod
16 use bMatrixHI_mod
17 use lamBmatrixHI_mod
18 use bMatrixEnsemble_mod
19 use bMatrixChem_mod
20 use bMatrixDiff_mod
21 use controlVector_mod
22 use verticalCoord_mod
23 use gridStateVector_mod
24 use horizontalCoord_mod
25 use utilities_mod
26 implicit none
27 save
28 private
29
30 ! public procedures
31 public :: bmat_setup, bmat_finalize, bmat_sqrtB, bmat_sqrtBT
32 public :: bmat_reduceToMPILocal, bmat_reduceToMPILocal_r4, bmat_expandToMPIGlobal, bmat_expandToMPIGlobal_r4
33
34 logical :: globalGrid = .true.
35
36 integer, parameter :: numMasterBmat = 4
37 character(len=4), parameter :: masterBmatTypeList (numMasterBmat) = (/'HI' , 'ENS' , 'CHM' , 'DIFF' /)
38 character(len=8), parameter :: masterBmatLabelList(numMasterBmat) = (/'B_HI', 'B_ENS', 'B_CHM', 'B_DIFF'/)
39 logical, parameter :: masterbmatIs3dList (numMasterBmat) = (/.true., .false., .true. , .true. /)
40
41 integer :: numBmat
42 integer, parameter :: numBmatMax = 50
43
44 character(len=4) :: bmatTypeList (numBmatMax)
45 character(len=9) :: bmatLabelList (numBmatMax)
46 integer :: bmatInstanceID(numBmatMax)
47 logical :: bmatIs3dList (numBmatMax)
48 logical :: bmatActive (numBmatMax)
49
50contains
51
52 !--------------------------------------------------------------------------
53 ! bmat_setup
54 !--------------------------------------------------------------------------
55 subroutine bmat_setup(hco_anl, hco_core, vco_anl)
56 !
57 !:Purpose: To initialize the analysis Background term for the specific
58 ! analysis configuration used.
59 !
60 implicit none
61
62 ! Arguments:
63 type(struct_vco), pointer, intent(in) :: vco_anl
64 type(struct_hco), pointer, intent(in) :: hco_anl
65 type(struct_hco), pointer, intent(in) :: hco_core
66
67 ! Locals:
68 integer, allocatable :: cvDimPerInstance(:)
69 integer :: cvdim
70 integer :: masterBmatIndex, bMatInstanceIndex, nBmatInstance, bmatIndex
71 character(len=2) :: bMatInstanceIndexString
72 character(len=3) :: bMatExtraLabel
73 logical :: active
74
75 call utl_tmg_start(50,'--Bmatrix')
76
77 !
78 !- 1. Get/Check the analysis grid info
79 !
80
81 !- 1.1 Horizontal Grid info
82 globalGrid = hco_anl%global
83
84 !
85 !- 2. Setup the B matrices
86 !
87 numBmat = 0
88
89 do masterBmatIndex = 1, numMasterBmat
90
91 select case( trim(masterBmatTypeList(masterBmatIndex)) )
92 case ('HI')
93
94 !- 2.1 Time-Mean Homogeneous and Isotropic...
95 nBmatInstance = 1 ! hardwired
96 allocate(cvdimPerInstance(nBmatInstance))
97
98 if ( globalGrid ) then
99 write(*,*)
100 write(*,*) 'Setting up the modular GLOBAL HI covariances...'
101 call bhi_Setup( hco_anl, vco_anl, & ! IN
102 cvdim ) ! OUT
103 else
104 write(*,*)
105 write(*,*) 'Setting up the modular LAM HI covariances...'
106 call lbhi_Setup( hco_anl, hco_core, vco_anl, & ! IN
107 cvdim ) ! OUT
108 end if
109
110 cvdimPerInstance(1) = cvdim
111
112 case ('ENS')
113
114 !- 2.2 Flow-dependent Ensemble-Based
115 write(*,*)
116 write(*,*) 'Setting up the modular ENSEMBLE covariances...'
117 call ben_Setup( hco_anl, hco_core, vco_anl, & ! IN
118 cvdimPerInstance ) ! OUT
119
120 nBmatInstance = size(cvdimPerInstance)
121
122 case ('CHM')
123
124 !- 2.3 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
125 nBmatInstance = 1 ! hardwired
126 allocate(cvdimPerInstance(nBmatInstance))
127
128 if ( globalGrid ) then
129 write(*,*)
130 write(*,*) 'Setting up the modular GLOBAL HI-chm covariances...'
131 call bchm_SetupCH( hco_anl, vco_anl, & ! IN
132 cvdim ) ! OUT
133 else
134 cvdim=0
135 end if
136
137 cvdimPerInstance(1) = cvdim
138
139 case ('DIFF')
140
141 !- 2.4 Covariances modelled using a diffusion operator.
142 nBmatInstance = 1 ! hardwired
143 allocate(cvdimPerInstance(nBmatInstance))
144
145 write(*,*)
146 write(*,*) 'Setting up the modular DIFFUSION covariances...'
147 call bdiff_Setup( hco_anl, vco_anl, & ! IN
148 cvdim ) ! OUT
149
150 cvdimPerInstance(1) = cvdim
151
152 case default
153
154 call utl_abort( 'bmat_setup: requested bmatrix type does not exist ' // trim(masterBmatTypeList(masterBmatIndex)) )
155
156 end select
157
158 !- 2.5 Append the info to the B matrix info arrays and setup the proper control sub-vectors
159 do bMatInstanceIndex = 1, nBmatInstance
160
161 numBmat = numBmat + 1
162
163 if (nBmatInstance == 1) then
164 bMatExtraLabel= ""
165 else
166 write(bMatInstanceIndexString,'(I2.2)') bMatInstanceIndex
167 bMatExtraLabel="_"//trim(bMatInstanceIndexString)
168 end if
169 bmatLabelList (numBmat) = trim(masterbmatLabelList(masterBmatIndex))//trim(bMatExtraLabel)
170
171 bmatTypeList (numBmat) = masterBmatTypeList(masterBmatIndex)
172 bmatIs3dList (numBmat) = masterbmatIs3dList(masterBmatIndex)
173 bmatInstanceID(numBmat) = bMatInstanceIndex
174
175 call cvm_setupSubVector(bmatLabelList(numBmat), bmatTypeList(numBmat), cvdimPerInstance(bMatInstanceIndex))
176
177 end do
178
179 deallocate(cvdimPerInstance)
180
181 end do
182
183 !
184 !- 3. Print a summary and set the active B matrices array
185 !
186 write(*,*)
187 write(*,*) " bmat_setup SUMMARY, number of B matrices found = ", numBmat
188 do bmatIndex = 1, numBmat
189 write(*,*) " B matrix #", bmatIndex
190 active = cvm_subVectorExists(bmatLabelList(bmatIndex))
191 if (active) then
192 write(*,*) " ACTIVE"
193 else
194 write(*,*) " NOT USED"
195 end if
196 write(*,*) " -> label = ", bmatLabelList (bmatIndex)
197 write(*,*) " -> type = ", bmatTypeList (bmatIndex)
198 if (active) then
199 write(*,*) " -> is 3D = ", bmatIs3dList (bmatIndex)
200 write(*,*) " -> instance ID = ", bmatInstanceID(bmatIndex)
201 end if
202 bmatActive(bmatIndex) = active
203 end do
204
205 !
206 !- 4. Check if the control vector has zero length (i.e. no B matrix is active)
207 !
208 if (cvm_nvadim_mpiglobal == 0) then
209 call utl_abort('bmat_setup: The control vector has zero length. No B matrix has been set up.')
210 end if
211
212 call utl_tmg_stop(50)
213
214 end subroutine bmat_setup
215
216 !--------------------------------------------------------------------------
217 ! bmat_sqrtB
218 !--------------------------------------------------------------------------
219 subroutine bmat_sqrtB(controlVector, cvdim, statevector, &
220 useFSOFcst_opt, stateVectorRef_opt)
221 !
222 !:Purpose: To transform model state from control-vector space to grid-point
223 ! space.
224 !
225 !
226 implicit none
227
228 ! Arguments:
229 integer , intent(in) :: cvdim
230 real(8) , intent(in) :: controlVector(cvdim)
231 type(struct_gsv) , intent(inout) :: statevector
232 logical, optional , intent(in) :: useFSOFcst_opt
233 type(struct_gsv), optional, intent(in) :: stateVectorRef_opt
234
235 ! Locals:
236 integer :: bmatIndex
237 real(8),pointer :: subVector(:)
238 type(struct_gsv) :: statevector_temp
239 character(len=4), pointer :: varNames(:)
240
241 call utl_tmg_start(50,'--Bmatrix')
242
243 !
244 !- 1. Set analysis increment to zero and allocate a temporary statevector
245 !
246 nullify(varNames)
247 call gsv_zero( statevector )
248 call gsv_varNamesList(varNames, statevector)
249 call gsv_allocate( statevector_temp, statevector%numStep, &
250 gsv_getHco(statevector), gsv_getVco(statevector), &
251 dataKind_opt=gsv_getDataKind(statevector), &
252 mpi_local_opt=.true., varNames_opt=varNames )
253 deallocate(varNames)
254
255 !
256 !- 2. Compute the analysis increment
257 !
258 bmat_loop: do bmatIndex = 1, numBmat
259
260 if ( .not. bmatActive(bmatIndex) ) cycle bmat_loop
261
262 subVector => cvm_getSubVector( controlVector, bmatLabelList(bmatIndex) )
263 call gsv_zero( statevector_temp )
264
265 select case( trim(bmatTypeList(bmatIndex)) )
266 case ('HI')
267
268 !- 2.1 Time-Mean Homogeneous and Isotropic...
269 call utl_tmg_start(52,'----B_HI_TL')
270 if ( globalGrid ) then
271 call bhi_bsqrt( subVector, & ! IN
272 statevector_temp, & ! OUT
273 stateVectorRef_opt ) ! IN
274 else
275 call lbhi_bSqrt( subVector, & ! IN
276 statevector_temp, & ! OUT
277 stateVectorRef_opt ) ! IN
278 end if
279 call utl_tmg_stop(52)
280
281 case ('CHM')
282
283 !- 2.2 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
284 call utl_tmg_start(68,'----B_CHM_TL')
285 if ( globalGrid ) then
286 call bchm_bsqrt( subVector, & ! IN
287 statevector_temp, & ! OUT
288 stateVectorRef_opt ) ! IN
289 end if
290 call utl_tmg_stop(68)
291
292 case ('DIFF')
293
294 !- 2.3 Covariances modelled using a diffusion operator.
295 call utl_tmg_start(66,'----B_DIFF_TL')
296 call bdiff_bsqrt( subVector, & ! IN
297 statevector_temp ) ! OUT
298 call utl_tmg_stop(66)
299
300 case ('ENS')
301
302 !- 2.4 Flow-dependent Ensemble-Based
303 call utl_tmg_start(57,'----B_ENS_TL')
304 call ben_bsqrt( bmatInstanceID(bmatIndex), subVector, & ! IN
305 statevector_temp, & ! OUT
306 useFSOFcst_opt, stateVectorRef_opt ) ! IN
307 call utl_tmg_stop(57)
308
309 end select
310
311 ! Make latest increment contribution 4D, if necessary
312 if ( bmatIs3dList(bmatIndex) ) call gsv_3dto4d( statevector_temp )
313
314 ! Add latest contribution to total increment in statevector
315 call gsv_add( statevector_temp, statevector )
316
317 end do bmat_loop
318
319 call gsv_deallocate( statevector_temp )
320
321 call utl_tmg_stop(50)
322
323 end subroutine bmat_sqrtB
324
325 !--------------------------------------------------------------------------
326 ! bmat_sqrtBT
327 !--------------------------------------------------------------------------
328 subroutine bmat_sqrtBT(controlVector, cvdim, statevector, &
329 useFSOFcst_opt, stateVectorRef_opt)
330 !
331 !:Purpose: To transform model state from grid-point space to
332 ! error-covariance space.
333 !
334 implicit none
335
336 ! Arguments:
337 integer , intent(in) :: cvdim
338 real(8) , intent(inout) :: controlVector(cvdim)
339 type(struct_gsv), intent(in) :: statevector
340 logical,optional, intent(in) :: useFSOFcst_opt
341 type(struct_gsv), optional, intent(in) :: stateVectorRef_opt
342
343 ! Locals:
344 integer :: bmatIndex
345 real(8),pointer :: subVector(:)
346 type(struct_gsv) :: statevector_temp
347 character(len=4), pointer :: varNames(:)
348
349 call utl_tmg_start(50,'--Bmatrix')
350
351 nullify(varNames)
352 call gsv_varNamesList(varNames, statevector)
353 call gsv_allocate( statevector_temp, statevector%numStep, &
354 gsv_getHco(statevector), gsv_getVco(statevector), &
355 dataKind_opt=gsv_getDataKind(statevector), &
356 mpi_local_opt=.true., varNames_opt=varNames )
357 deallocate(varNames)
358
359 ! Process components in opposite order as forward calculation
360 bmat_loop: do bmatIndex = numBmat, 1, -1
361
362 if ( .not. bmatActive(bmatIndex) ) cycle bmat_loop
363
364 subVector => cvm_getSubVector( controlVector, bmatLabelList(bmatIndex) )
365 subVector(:) = 0.0d0
366
367 ! Adjoint of converting statevector from 3D to 4D
368 call gsv_copy( statevector, statevector_temp )
369 if ( bmatIs3dList(bmatIndex) ) call gsv_3dto4dAdj( statevector_temp )
370
371 select case( trim(bmatTypeList(bmatIndex)) )
372 case ('ENS')
373
374 !- 2.1 Flow-dependent Ensemble-Based
375 call utl_tmg_start(61,'----B_ENS_AD')
376
377 call ben_bsqrtad( bmatInstanceID(bmatIndex), statevector_temp, & ! IN
378 subVector, & ! OUT
379 useFSOFcst_opt, stateVectorRef_opt ) ! IN
380 call utl_tmg_stop(61)
381
382 case ('DIFF')
383
384 !- 2.2 Covariances modelled using a diffusion operator.
385 call utl_tmg_start(67,'----B_DIFF_AD')
386 call bdiff_bsqrtad( statevector_temp, & ! IN
387 subVector ) ! OUT
388 call utl_tmg_stop(67)
389
390 case ('CHM')
391
392 !- 2.3 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
393 call utl_tmg_start(69,'----B_CHM_AD')
394 if ( globalGrid ) then
395 call bchm_bsqrtad( statevector_temp, & ! IN
396 subVector, & ! OUT
397 stateVectorRef_opt ) ! IN
398 end if
399 call utl_tmg_stop(69)
400
401 case ('HI')
402
403 !- 2.4 Time-Mean Homogeneous and Isotropic...
404 call utl_tmg_start(53,'----B_HI_AD')
405 if ( globalGrid ) then
406 call bhi_bsqrtad( statevector_temp, & ! IN
407 subVector, & ! OUT
408 stateVectorRef_opt ) ! IN
409 else
410 call lbhi_bSqrtAdj( statevector_temp, & ! IN
411 subVector, & ! OUT
412 stateVectorRef_opt ) ! IN
413 end if
414 call utl_tmg_stop(53)
415
416 end select
417
418 end do bmat_loop
419
420 call gsv_deallocate( statevector_temp )
421
422 call utl_tmg_stop(50)
423
424 end subroutine bmat_sqrtBT
425
426 !--------------------------------------------------------------------------
427 ! bmat_finalize
428 !--------------------------------------------------------------------------
429 subroutine bmat_finalize()
430 !
431 !:Purpose: To release memory used by B matrices.
432 !
433 implicit none
434
435 call bhi_finalize()
436 call ben_finalize()
437 call bchm_finalize()
438 call lbhi_finalize()
439 call bdiff_finalize()
440
441 end subroutine bmat_finalize
442
443 !--------------------------------------------------------------------------
444 ! bmat_reduceToMPILocal
445 !--------------------------------------------------------------------------
446 subroutine bmat_reduceToMPILocal(cv_mpilocal,cv_mpiglobal)
447 !
448 !:Purpose: To distribute MPI_global control vector from task 0 to all tasks
449 ! where the arguments are real(8)'s
450 !
451 implicit none
452
453 ! Arguments:
454 real(8), intent(inout) :: cv_mpilocal(:)
455 real(8), intent(in) :: cv_mpiglobal(:)
456
457 ! Locals:
458 integer :: bmatIndex
459 real(8), pointer :: subVector_mpilocal(:), subVector_mpiglobal(:)
460 real(8), target :: dummyVector(1)
461
462 bmat_loop: do bmatIndex = 1, numBmat
463
464 if ( .not. bmatActive(bmatIndex) ) cycle bmat_loop
465
466 subVector_mpilocal => cvm_getSubVector( cv_mpilocal, bmatLabelList(bmatIndex) )
467 if ( mmpi_myid == 0 ) then
468 subVector_mpiglobal => cvm_getSubVector_mpiglobal( cv_mpiglobal, bmatLabelList(bmatIndex) )
469 else
470 subVector_mpiglobal => dummyVector
471 end if
472
473 select case( trim(bmatTypeList(bmatIndex)) )
474 case ('HI')
475
476 !- 2.1 Time-Mean Homogeneous and Isotropic...
477 if ( globalGrid ) then
478 call bhi_reduceToMPILocal( subVector_mpilocal,subVector_mpiglobal )
479 else
480 call lbhi_reduceToMPILocal( subVector_mpilocal,subVector_mpiglobal )
481 end if
482
483 case ('CHM')
484
485 !- 2.2 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
486 if ( globalGrid ) then
487 call bchm_reduceToMPILocal( subVector_mpilocal, subVector_mpiglobal )
488 end if
489
490 case ('DIFF')
491
492 !- 2.3 Covariances modelled using a diffusion operator.
493 call bdiff_reduceToMPILocal( subVector_mpilocal, subVector_mpiglobal )
494
495 case ('ENS')
496
497 !- 2.4 Flow-dependent Ensemble-Based
498 call ben_reduceToMPILocal(subVector_mpilocal, subVector_mpiglobal, bmatInstanceID(bmatIndex))
499
500 end select
501
502 end do bmat_loop
503
504 end subroutine bmat_reduceToMPILocal
505
506 !--------------------------------------------------------------------------
507 ! bmat_reduceToMPILocal_r4
508 !--------------------------------------------------------------------------
509 subroutine bmat_reduceToMPILocal_r4(cv_mpilocal,cv_mpiglobal)
510 !
511 !:Purpose: To distribute MPI_global control vector from task 0 to all tasks
512 ! where the arguments are real(4)'s.
513 !
514 implicit none
515
516 ! Arguments:
517 real(4), intent(inout) :: cv_mpilocal(:)
518 real(4), intent(in) :: cv_mpiglobal(:)
519
520 ! Locals:
521 integer :: bmatIndex
522 real(4), pointer :: subVector_mpilocal(:), subVector_mpiglobal(:)
523 real(4), target :: dummyVector_r4(1)
524
525 bmat_loop: do bmatIndex = 1, numBmat
526
527 if ( .not. bmatActive(bmatIndex) ) cycle bmat_loop
528
529 subVector_mpilocal => cvm_getSubVector_r4( cv_mpilocal, bmatLabelList(bmatIndex) )
530 if ( mmpi_myid == 0 ) then
531 subVector_mpiglobal => cvm_getSubVector_mpiglobal_r4( cv_mpiglobal, bmatLabelList(bmatIndex) )
532 else
533 subVector_mpiglobal => dummyVector_r4
534 end if
535
536 select case( trim(bmatTypeList(bmatIndex)) )
537 case ('HI')
538
539 !- 2.1 Time-Mean Homogeneous and Isotropic...
540 if ( globalGrid ) then
541 call bhi_reduceToMPILocal_r4( subVector_mpilocal,subVector_mpiglobal )
542 else
543 call lbhi_reduceToMPILocal_r4( subVector_mpilocal,subVector_mpiglobal )
544 end if
545
546 case ('CHM')
547
548 !- 2.2 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
549 if ( globalGrid ) then
550 call bchm_reduceToMPILocal_r4( subVector_mpilocal, subVector_mpiglobal )
551 end if
552
553 case ('DIFF')
554
555 !- 2.3 Covariances modelled using a diffusion operator.
556 call utl_abort('bmat_reduceToMPILocal_r4: not yet implemented for bMatrixDiff')
557 !call bdiff_reduceToMPILocal_r4( subVector_mpilocal, subVector_mpiglobal )
558
559 case ('ENS')
560
561 !- 2.4 Flow-dependent Ensemble-Based
562 call ben_reduceToMPILocal_r4(subVector_mpilocal, subVector_mpiglobal, bmatInstanceID(bmatIndex))
563
564 end select
565
566 end do bmat_loop
567
568 end subroutine bmat_reduceToMPILocal_r4
569
570 !--------------------------------------------------------------------------
571 ! bmat_expandToMPIGlobal
572 !--------------------------------------------------------------------------
573 subroutine bmat_expandToMPIGlobal(cv_mpilocal,cv_mpiglobal)
574 !
575 !:Purpose: To gather control vector from all tasks to task 0 where the
576 ! arguments are real(8)'s.
577 !
578 implicit none
579
580 ! Arguments:
581 real(8), intent(in) :: cv_mpilocal(:)
582 real(8), intent(inout) :: cv_mpiglobal(:)
583
584 ! Locals:
585 integer :: bmatIndex
586 real(8), pointer :: subVector_mpilocal(:), subVector_mpiglobal(:)
587 real(8), target :: dummyVector(1)
588
589 bmat_loop: do bmatIndex = 1, numBmat
590
591 if ( .not. bmatActive(bmatIndex) ) cycle bmat_loop
592
593 subVector_mpilocal => cvm_getSubVector( cv_mpilocal, bmatLabelList(bmatIndex) )
594 if ( mmpi_myid == 0 ) then
595 subVector_mpiglobal => cvm_getSubVector_mpiglobal( cv_mpiglobal, bmatLabelList(bmatIndex) )
596 else
597 subVector_mpiglobal => dummyVector
598 end if
599
600 select case( trim(bmatTypeList(bmatIndex)) )
601 case ('HI')
602
603 !- 2.1 Time-Mean Homogeneous and Isotropic...
604 if ( globalGrid ) then
605 call bhi_expandToMPIGlobal( subVector_mpilocal,subVector_mpiglobal )
606 else
607 call lbhi_expandToMPIGlobal( subVector_mpilocal,subVector_mpiglobal )
608 end if
609
610 case ('CHM')
611
612 !- 2.2 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
613 if ( globalGrid ) then
614 call bchm_expandToMPIGlobal( subVector_mpilocal, subVector_mpiglobal )
615 end if
616
617 case ('DIFF')
618
619 !- 2.3 Covariances modelled using a diffusion operator.
620 !call bdiff_expandToMPIGlobal( subVector_mpilocal, subVector_mpiglobal )
621
622 case ('ENS')
623
624 !- 2.4 Flow-dependent Ensemble-Based
625 call ben_expandToMPIGlobal(subVector_mpilocal, subVector_mpiglobal, bmatInstanceID(bmatIndex) )
626
627 end select
628
629 end do bmat_loop
630
631 end subroutine bmat_expandToMPIGlobal
632
633 !--------------------------------------------------------------------------
634 ! bmat_expandToMPIGlobal_r4
635 !--------------------------------------------------------------------------
636 subroutine bmat_expandToMPIGlobal_r4(cv_mpilocal,cv_mpiglobal)
637 !
638 !:Purpose: To gather control vector from all tasks to task 0 where the
639 ! arguments are real(4)'s.
640 !
641 implicit none
642
643 ! Arguments:
644 real(4), intent(in) :: cv_mpilocal(:)
645 real(4), intent(inout) :: cv_mpiglobal(:)
646
647 ! Locals:
648 integer :: bmatIndex
649 real(4), pointer :: subVector_mpilocal(:), subVector_mpiglobal(:)
650 real(4), target :: dummyVector_r4(1)
651
652 bmat_loop: do bmatIndex = 1, numBmat
653
654 if ( .not. bmatActive(bmatIndex) ) cycle bmat_loop
655
656 subVector_mpilocal => cvm_getSubVector_r4( cv_mpilocal, bmatLabelList(bmatIndex) )
657 if ( mmpi_myid == 0 ) then
658 subVector_mpiglobal => cvm_getSubVector_mpiglobal_r4( cv_mpiglobal, bmatLabelList(bmatIndex) )
659 else
660 subVector_mpiglobal => dummyVector_r4
661 end if
662
663 select case( trim(bmatTypeList(bmatIndex)) )
664 case ('HI')
665
666 !- 2.1 Time-Mean Homogeneous and Isotropic...
667 if ( globalGrid ) then
668 call bhi_expandToMPIGlobal_r4( subVector_mpilocal,subVector_mpiglobal )
669 else
670 call lbhi_expandToMPIGlobal_r4( subVector_mpilocal,subVector_mpiglobal )
671 end if
672
673 case ('CHM')
674
675 !- 2.2 Static (Time-Mean Homogeneous and Isotropic) covariances for constituents
676 if ( globalGrid ) then
677 call bchm_expandToMPIGlobal_r4( subVector_mpilocal, subVector_mpiglobal )
678 end if
679
680 case ('DIFF')
681
682 !- 2.3 Covariances modelled using a diffusion operator.
683 !call bdiff_expandToMPIGlobal_r4( subVector_mpilocal, subVector_mpiglobal )
684
685 case ('ENS')
686
687 !- 2.4 Flow-dependent Ensemble-Based
688 call ben_expandToMPIGlobal_r4(subVector_mpilocal, subVector_mpiglobal, bmatInstanceID(bmatIndex) )
689
690 end select
691
692 end do bmat_loop
693
694 end subroutine bmat_expandToMPIGlobal_r4
695
696end module bMatrix_mod