1module columnData_mod
2 ! MODULE columnData_mod (prefix='col' category='6. High-level data objects')
3 !
4 !:Purpose: A derived type and related procedures for storing and manipulating
5 ! vertical columns of analysis variables on model or analysis grid
6 ! levels. These columns are generally produced by horizontally
7 ! interpolating a gridStateVector object to the observation
8 ! locations.
9 !
10 use midasMpi_mod
11 use varNameList_mod
12 use verticalCoord_mod
13 use mathPhysConstants_mod
14 use utilities_mod
15
16 implicit none
17 save
18 private
19
20 ! public variables and types
21 public :: col_rhumin, col_minValVarKindCH, struct_columnData
22
23 ! public subroutines and functions
24 public :: col_setup, col_allocate, col_deallocate
25 public :: col_varExist, col_getOffsetFromVarno
26 public :: col_getNumLev, col_getNumCol, col_getVarNameFromK
27 public :: col_getPressure, col_getHeight, col_setHeightSfc
28 public :: col_zero, col_getAllColumns, col_getColumn, col_getElem, col_getVco, col_setVco
29 public :: col_getLevIndexFromVarLevIndex
30
31 type struct_columnData
32 integer :: nk, numCol
33 logical :: allocated=.false.
34 logical :: mpi_local
35 real(8), pointer :: all(:,:)
36 real(8), pointer :: heightSfc(:)
37 real(8), pointer :: oltv(:,:,:) ! Tangent linear operator of virtual temperature
38 integer, pointer :: varOffset(:),varNumLev(:)
39 logical :: varExistList(vnl_numVarMax)
40 type(struct_vco), pointer :: vco => null()
41 logical :: addHeightSfcOffset = .false.
42 real(8), pointer :: lat(:)
43 end type struct_columnData
44
45 real(8) :: col_rhumin
46 logical :: varExistList(vnl_numvarmax)
47
48 ! Minimum values for variables of CH kind
49 real(8) :: col_minValVarKindCH(vnl_numVarMax)
50
51 ! Namelist variables
52 real(8) :: rhumin ! minimum humidity value imposed after interpolation to columns
53 logical :: addHeightSfcOffset ! choose to add non-zero height offset to diagnostic (sfc) levels
54 real(8) :: minValVarKindCH(vnl_numVarMax) ! variable-dependent minimum value applied to chemistry variables
55
56contains
57
58 !--------------------------------------------------------------------------
59 ! col_setup
60 !--------------------------------------------------------------------------
61 subroutine col_setup
62 implicit none
63
64 ! Locals:
65 integer :: varIndex, loopIndex
66 integer :: fnom,fclos,nulnam,ierr
67 integer :: numVar3D, numVar2D, numVarOther
68
69 ! Namelist variables (local)
70 character(len=4) :: anlvar(vnl_numvarmax) ! list of state variable names
71 character(len=8) :: anltime_bin ! can be 'MIDDLE', 'FIRST' or 'LAST'
72 logical :: conversionVarKindCHtoMicrograms ! activate unit conversion for CH variables
73 logical :: abortOnMpiImbalance ! choose to abort program when MPI imbalance is too large
74
75 namelist /namstate/anlvar,rhumin,anltime_bin,addHeightSfcOffset,conversionVarKindCHtoMicrograms, &
76 minValVarKindCH, abortOnMpiImbalance
77
78 if(mmpi_myid == 0) write(*,*) 'col_setup: List of known (valid) variable names'
79 if(mmpi_myid == 0) write(*,*) 'col_setup: varNameList3D=',vnl_varNameList3D
80 if(mmpi_myid == 0) write(*,*) 'col_setup: varNameList2D=',vnl_varNameList2D
81 if(mmpi_myid == 0) write(*,*) 'col_setup: varNameList =',vnl_varNameList
82
83 ! Read NAMELIST NAMSTATE to find which fields are needed
84
85 anlvar(:) = ' '
86 rhumin = MPC_MINIMUM_HU_R8
87 anltime_bin = 'MIDDLE'
88 addHeightSfcOffset = .false.
89 conversionVarKindCHtoMicrograms = .false.
90 minValVarKindCH(:) = mpc_missingValue_r8
91 abortOnMpiImbalance = .true.
92
93 nulnam=0
94 ierr=fnom(nulnam,'./flnml','FTN+SEQ+R/O',0)
95 read(nulnam,nml=namstate,iostat=ierr)
96 if(ierr.ne.0) call utl_abort('col_setup: Error reading namelist')
97 if(mmpi_myid == 0) write(*,nml=namstate)
98 ierr=fclos(nulnam)
99
100 col_rhumin = rhumin
101 col_minValVarKindCH(:)=minValVarKindCH(:)
102
103 if( varneed('Z_T') .or. varneed('Z_M') ) then
104 call utl_abort('col_setup: height can not be specified as analysis variable in namelist!')
105 end if
106 if( varneed('P_T') .or. varneed('P_M') ) then
107 call utl_abort('col_setup: pressure can not be specified as analysis variable in namelist!')
108 end if
109
110 numVar3D = 0
111 numVar2D = 0
112 numVarOther = 0
113
114 do varIndex = 1, vnl_numvarmax3D
115 if (varneed(vnl_varNameList3D(varIndex))) then
116 varExistList(varIndex) = .true.
117 numVar3D = numVar3D + 1
118 else
119 varExistList(varIndex) = .false.
120 end if
121 end do
122
123 do varIndex = 1, vnl_numvarmax2D
124 if (varneed(vnl_varNameList2D(varIndex))) then
125 varExistList(varIndex+vnl_numvarmax3D) = .true.
126 numVar2D = numVar2D + 1
127 else
128 varExistList(varIndex+vnl_numvarmax3D) = .false.
129 end if
130 end do
131
132 do varIndex = 1, vnl_numvarmaxOther
133 if (varneed(vnl_varNameListOther(varIndex))) then
134 varExistList(varIndex+vnl_numvarmax3D+vnl_numvarmax2D) = .true.
135 numVarOther = numVarOther + 1
136 else
137 varExistList(varIndex+vnl_numvarmax3D+vnl_numvarmax2D) = .false.
138 end if
139 end do
140
141 ! Setup to assign min values to apply
142
143 ! Check for input values only for variables of CH kind
144 do varIndex = 1, vnl_numvarmax
145 if ( trim(AnlVar(varIndex)) == '' ) exit
146 if ( vnl_varKindFromVarname(AnlVar(varIndex)) == 'CH' ) then
147 if ( minValVarKindCH(varIndex) < 0.99d0 * MPC_missingValue_R8 ) then
148 if ( trim(AnlVar(varIndex)) == 'AF' .or. trim(AnlVar(varIndex)) == 'AC' ) then
149 ! Set for particulate matter in micrograms/cm^3
150 minValVarKindCH(varIndex) = MPC_MINIMUM_PM_R8
151 else
152 ! Set for concentrations in micrograms/kg
153 minValVarKindCH(varIndex) = MPC_MINIMUM_CH_R8
154 end if
155 end if
156 end if
157 end do
158
159 ! Assign min values to apply
160 col_minValVarKindCH(:) = MPC_missingValue_R8
161 do varIndex = 1, vnl_numvarmax
162 if ( varExistList(varIndex) ) then
163 do loopIndex = 1, vnl_numvarmax
164 if ( trim(AnlVar(loopIndex)) == '' ) exit
165 if ( trim(vnl_varNameList(varIndex)) == trim(AnlVar(loopIndex)) ) &
166 col_minValVarKindCH(varIndex) = minValVarKindCH(loopIndex)
167 end do
168 end if
169 end do
170
171 if(mmpi_myid == 0) write(*,*) 'col_setup: numVar3D (no Z_T/Z_M/P_T/P_M included), numVar2D, numVarOther = ', numVar3D, numVar2D, numVarOther
172 if(mmpi_myid == 0) write(*,*) 'col_setup: varExistList (no Z_T/Z_M/P_T/P_M included) = ',varExistList
173
174 contains
175
176 logical function varneed(varName)
177 character(len=*) :: varName
178 integer :: jvar
179
180 varneed = .false.
181 do jvar = 1, vnl_numVarMax
182 if (trim(varName) == trim(anlvar(jvar))) then
183 varneed = .true.
184 end if
185 end do
186
187 end function varneed
188
189 end subroutine col_setup
190
191 !--------------------------------------------------------------------------
192 ! col_zero
193 !--------------------------------------------------------------------------
194 subroutine col_zero(column)
195 implicit none
196
197 ! Arguments:
198 type(struct_columnData), intent(inout) :: column
199
200 if (column%numCol > 0) then
201 column%all(:,:) = 0.0d0
202 column%heightSfc(:) = 0.0d0
203 end if
204
205 end subroutine col_zero
206
207 !--------------------------------------------------------------------------
208 ! col_allocate
209 !--------------------------------------------------------------------------
210 subroutine col_allocate(column, numCol, mpiLocal_opt, beSilent_opt, setToZero_opt, varNames_opt)
211 implicit none
212
213 ! Arguments:
214 type(struct_columnData), intent(inout) :: column
215 integer, intent(in) :: numCol
216 logical, optional, intent(in) :: mpiLocal_opt
217 logical, optional, intent(in) :: beSilent_opt
218 logical, optional, intent(in) :: setToZero_opt
219 character(len=*), optional, intent(in) :: varNames_opt(:)
220
221 ! Locals:
222 integer :: iloc, varIndex, varIndex2, numVar
223 logical :: beSilent, setToZero
224
225 if ( present(beSilent_opt) ) then
226 beSilent = beSilent_opt
227 else
228 beSilent = .false.
229 end if
230
231 if ( present(setToZero_opt) ) then
232 setToZero = setToZero_opt
233 else
234 setToZero = .true.
235 end if
236
237 if ( present(varNames_opt) ) then
238 column%varExistList(:) = .false.
239 numVar = size( varNames_opt )
240 do varIndex2 = 1, numVar
241 varIndex = vnl_varListIndex(varNames_opt(varIndex2))
242 column%varExistList(varIndex) = .true.
243 end do
244 else
245 ! set the variable list using the global ExistList
246 column%varExistList(:) = varExistList(:)
247 end if
248
249 if ( column%varExistList(vnl_varListIndex('TT')) .and. &
250 column%varExistList(vnl_varListIndex('HU')) .and. &
251 column%varExistList(vnl_varListIndex('P0')) ) then
252 if ( col_getNumLev(column,'TH') > 0 ) column%varExistList(vnl_varListIndex('Z_T')) = .true.
253 if ( col_getNumLev(column,'MM') > 0 ) column%varExistList(vnl_varListIndex('Z_M')) = .true.
254 end if
255
256 if ( column%varExistList(vnl_varListIndex('P0')) ) then
257 if ( col_getNumLev(column,'TH') > 0 ) column%varExistList(vnl_varListIndex('P_T')) = .true.
258 if ( col_getNumLev(column,'MM') > 0 ) column%varExistList(vnl_varListIndex('P_M')) = .true.
259 end if
260
261 ! add P0LS to the varExistList if vcode=5100
262 if (column%vco%vcode == 5100) then
263 column%varExistList(vnl_varListIndex('P0LS')) = .true.
264 end if
265
266 column%numCol = numCol
267 if(present(mpiLocal_opt)) then
268 column%mpi_local=mpiLocal_opt
269 else
270 column%mpi_local=.true.
271 if (mmpi_myid == 0 .and. .not.beSilent) write(*,*) 'col_allocate: assuming columnData is mpi-local'
272 end if
273
274 if(.not.column%vco%initialized) then
275 call utl_abort('col_allocate: VerticalCoord has not been initialized!')
276 end if
277
278 allocate(column%varOffset(vnl_numvarmax))
279 column%varOffset(:)=0
280 allocate(column%varNumLev(vnl_numvarmax))
281 column%varNumLev(:)=0
282
283 iloc = 0
284 do varIndex = 1, vnl_numvarmax3d
285 if(column%varExistList(varIndex)) then
286 column%varOffset(varIndex) = iloc
287 column%varNumLev(varIndex) = col_getNumLev(column,vnl_varLevelFromVarname(vnl_varNameList(varIndex)))
288 iloc = iloc + column%varNumLev(varIndex)
289 end if
290 end do
291 do varIndex2 = 1, vnl_numvarmax2d
292 varIndex = varIndex2+vnl_numvarmax3d
293 if(column%varExistList(varIndex)) then
294 column%varOffset(varIndex) = iloc
295 column%varNumLev(varIndex) = 1
296 iloc = iloc + 1
297 end if
298 end do
299 do varIndex2 = 1, vnl_numvarmaxOther
300 varIndex = varIndex2+vnl_numvarmax3d+vnl_numvarmax2d
301 if(column%varExistList(varIndex)) then
302 column%varOffset(varIndex) = iloc
303 column%varNumLev(varIndex) = col_getNumLev(column,'OT',vnl_varNameListOther(varIndex2))
304 iloc = iloc + column%varNumLev(varIndex)
305 end if
306 end do
307
308 if (iloc == 0) then
309 call utl_abort('col_allocate: Nothing to allocate')
310 end if
311
312 column%nk = iloc
313
314 if(column%numCol.le.0) then
315 if ( .not.beSilent ) write(*,*) 'col_allocate: number of columns is zero, not allocated'
316 else
317 allocate(column%all(column%nk,column%numCol))
318 if ( setToZero ) column%all(:,:)=0.0d0
319
320 allocate(column%heightSfc(column%numCol))
321 column%heightSfc(:)=0.0d0
322
323 allocate(column%oltv(2,col_getNumLev(column,'TH'),numCol))
324 if ( setToZero ) column%oltv(:,:,:)=0.0d0
325
326 allocate(column%lat(numCol))
327 if ( setToZero ) column%lat(:)=0.0d0
328 end if
329
330 if(mmpi_myid == 0 .and. .not.beSilent) write(*,*) 'col_allocate: column%nk = ', column%nk
331 if(mmpi_myid == 0 .and. .not.beSilent) write(*,*) 'col_allocate: varOffset=',column%varOffset
332 if(mmpi_myid == 0 .and. .not.beSilent) write(*,*) 'col_allocate: varNumLev=',column%varNumLev
333
334 column%addHeightSfcOffset = addHeightSfcOffset
335
336 column%allocated=.true.
337
338 end subroutine col_allocate
339
340 !--------------------------------------------------------------------------
341 ! col_deallocate
342 !--------------------------------------------------------------------------
343 subroutine col_deallocate(column)
344 implicit none
345
346 ! Arguments:
347 type(struct_columnData), intent(inout) :: column
348
349 deallocate(column%varOffset)
350 deallocate(column%varNumLev)
351
352 if(column%numCol.gt.0) then
353 deallocate(column%all)
354 deallocate(column%heightSfc)
355 deallocate(column%oltv)
356 end if
357
358 column%allocated=.false.
359
360 end subroutine col_deallocate
361
362 !--------------------------------------------------------------------------
363 ! col_varExist
364 !--------------------------------------------------------------------------
365 recursive function col_varExist(column_opt,varName) result(varExist)
366 implicit none
367
368 ! Arguments:
369 type(struct_columnData), optional, intent(in) :: column_opt
370 character(len=*), intent(in) :: varName
371 ! Result:
372 logical :: varExist
373
374 if (varName == 'Z_*') then
375 varExist = col_varExist(column_opt, 'Z_T') .and. &
376 col_varExist(column_opt, 'Z_M')
377 else if (varName == 'P_*') then
378 varExist = col_varExist(column_opt, 'P_T') .and. &
379 col_varExist(column_opt, 'P_M')
380 else
381 if ( present(column_opt) ) then
382 if ( column_opt%varExistList(vnl_varListIndex(varName)) ) then
383 varExist = .true.
384 else
385 varExist = .false.
386 end if
387 else
388 if ( varExistList(vnl_varListIndex(varName)) ) then
389 varExist = .true.
390 else
391 varExist = .false.
392 end if
393 end if
394
395 if (present(column_opt)) then
396 varExist = column_opt % varExistList(vnl_varListIndex(varName))
397 else
398 varExist = varExistList(vnl_varListIndex(varName))
399 end if
400 end if
401
402 end function col_varExist
403
404 !--------------------------------------------------------------------------
405 ! col_getOffsetFromVarno
406 !--------------------------------------------------------------------------
407 function col_getOffsetFromVarno(column,varnum,varNumberChm_opt,modelName_opt) result(offset)
408 implicit none
409
410 ! Arguments:
411 type(struct_columnData), intent(in) :: column
412 integer, intent(in) :: varnum
413 integer, optional, intent(in) :: varNumberChm_opt
414 character(len=*), optional, intent(in) :: modelName_opt
415 ! Result:
416 integer :: offset
417
418 offset=column%varOffset(vnl_varListIndex(vnl_varnameFromVarnum(varnum,varNumberChm_opt=varNumberChm_opt,modelName_opt=modelName_opt)))
419
420 end function col_getOffsetFromVarno
421
422 !--------------------------------------------------------------------------
423 ! col_getLevIndexFromVarLevIndex
424 !--------------------------------------------------------------------------
425 function col_getLevIndexFromVarLevIndex(column, varLevIndex) result(levIndex)
426 implicit none
427
428 ! Arguments:
429 type(struct_columnData), intent(in) :: column
430 integer, intent(in) :: varLevIndex
431 ! Result:
432 integer :: varIndex
433
434 ! Locals:
435 integer :: levIndex
436
437 do varIndex = 1, vnl_numvarmax
438 if ( column%varExistList(varIndex) ) then
439 if ( (varLevIndex >= (column%varOffset(varIndex) + 1)) .and. &
440 (varLevIndex <= (column%varOffset(varIndex) + column%varNumLev(varIndex))) ) then
441 levIndex = varLevIndex - column%varOffset(varIndex)
442 return
443 end if
444 end if
445 end do
446
447 write(*,*) 'col_getLevIndexFromVarLevIndex: varLevIndex out of range: ', varLevIndex
448 call utl_abort('col_getLevIndexFromVarLevIndex')
449
450 end function col_getLevIndexFromVarLevIndex
451
452 !--------------------------------------------------------------------------
453 ! col_getVarNameFromK
454 !--------------------------------------------------------------------------
455 function col_getVarNameFromK(column,kIndex) result(varName)
456 implicit none
457
458 ! Arguments:
459 type(struct_columnData), intent(in) :: column
460 integer, intent(in) :: kIndex
461 ! Result:
462 character(len=4) :: varName
463
464 ! Locals:
465 integer :: varIndex
466
467 do varIndex = 1, vnl_numvarmax
468 if ( column%varExistList(varIndex) ) then
469 if ( (kIndex >= (column%varOffset(varIndex) + 1)) .and. &
470 (kIndex <= (column%varOffset(varIndex) + column%varNumLev(varIndex))) ) then
471 varName = vnl_varNameList(varIndex)
472 return
473 end if
474 end if
475 end do
476
477 write(*,*) 'col_getVarNameFromK: kIndex out of range: ', kIndex
478 call utl_abort('col_getVarNameFromK')
479
480 end function col_getVarNameFromK
481
482 !--------------------------------------------------------------------------
483 ! col_getPressure
484 !--------------------------------------------------------------------------
485 function col_getPressure(column,ilev,headerIndex,varLevel) result(pressure)
486 implicit none
487
488 ! Arguments:
489 type(struct_columnData), intent(in) :: column
490 integer, intent(in) :: ilev
491 integer, intent(in) :: headerIndex
492 character(len=*), intent(in) :: varLevel
493 ! Result:
494 real(8) :: pressure
495
496 ! Locals:
497 integer :: ilev1
498
499 if (varLevel == 'TH' .and. col_varExist(column,'P_T')) then
500 ilev1 = 1 + column%varOffset(vnl_varListIndex('P_T'))
501 pressure = column%all(ilev1+ilev-1,headerIndex)
502 elseif (varLevel == 'MM' .and. col_varExist(column,'P_M') ) then
503 ilev1 = 1 + column%varOffset(vnl_varListIndex('P_M'))
504 pressure = column%all(ilev1+ilev-1,headerIndex)
505 else
506 call utl_abort('col_getPressure: Unknown variable type: ' // varLevel)
507 end if
508
509 end function col_getPressure
510
511 !--------------------------------------------------------------------------
512 ! col_getHeight
513 !--------------------------------------------------------------------------
514 function col_getHeight(column,ilev,headerIndex,varLevel) result(height)
515 implicit none
516
517 ! Arguments:
518 type(struct_columnData), intent(in) :: column
519 integer, intent(in) :: ilev
520 integer, intent(in) :: headerIndex
521 character(len=*), intent(in) :: varLevel
522 ! Result:
523 real(8) :: height
524
525 ! Locals:
526 integer :: ilev1
527
528 if (varLevel == 'TH') then
529 if (.not. col_varExist(column,'Z_T') ) then
530 call utl_abort('col_getHeight: Z_T not found!')
531 end if
532 ilev1 = 1 + column%varOffset(vnl_varListIndex('Z_T'))
533 height = column%all(ilev1+ilev-1,headerIndex)
534 else if (varLevel == 'MM') then
535 if (.not. col_varExist(column,'Z_M') ) then
536 call utl_abort('col_getHeight: Z_M not found!')
537 end if
538 ilev1 = 1 + column%varOffset(vnl_varListIndex('Z_M'))
539 height = column%all(ilev1+ilev-1,headerIndex)
540 else if (varLevel == 'SF' ) then
541 height = column%heightSfc(headerIndex)
542 else
543 call utl_abort('col_getHeight: unknown varLevel! ' // varLevel)
544 end if
545
546 end function col_getHeight
547
548 !--------------------------------------------------------------------------
549 ! col_setHeightsSfc
550 !--------------------------------------------------------------------------
551 subroutine col_setHeightSfc(column,headerIndex,height)
552 implicit none
553
554 ! Arguments:
555 type(struct_columnData), intent(inout) :: column
556 integer, intent(in) :: headerIndex
557 real(8), intent(in) :: height
558
559 column%heightSfc(headerIndex) = height
560
561 end subroutine col_setHeightSfc
562
563 !--------------------------------------------------------------------------
564 ! col_getAllColumns
565 !--------------------------------------------------------------------------
566 function col_getAllColumns(column,varName_opt) result(allColumns)
567 implicit none
568
569 ! Arguments:
570 type(struct_columnData), intent(in) :: column
571 character(len=*), optional, intent(in) :: varName_opt
572 ! Result:
573 real(8), pointer :: allColumns(:,:)
574
575 ! Locals:
576 integer :: ilev1,ilev2
577
578 if ( column%numCol > 0 ) then
579 if(present(varName_opt)) then
580 if ( col_varExist(column,varName_opt) ) then
581 ilev1 = column%varOffset(vnl_varListIndex(varName_opt))+1
582 ilev2 = ilev1 - 1 + column%varNumLev(vnl_varListIndex(varName_opt))
583 allColumns => column%all(ilev1:ilev2,:)
584 else
585 call utl_abort('col_getAllColumns: Unknown variable name! ' // varName_opt)
586 end if
587 else
588 allColumns => column%all(:,:)
589 end if
590 else
591 allColumns => null()
592 end if
593
594 end function col_getAllColumns
595
596 !--------------------------------------------------------------------------
597 ! col_getColumn
598 !--------------------------------------------------------------------------
599 function col_getColumn(column,headerIndex,varName_opt) result(onecolumn)
600 implicit none
601
602 ! Arguments:
603 type(struct_columnData), intent(in) :: column
604 integer, intent(in) :: headerIndex
605 character(len=*), optional, intent(in) :: varName_opt
606 ! Result:
607 real(8), pointer :: onecolumn(:)
608
609 ! Locals:
610 integer :: ilev1,ilev2
611
612 if(present(varName_opt)) then
613 if(col_varExist(column,varName_opt)) then
614 ilev1 = column%varOffset(vnl_varListIndex(varName_opt))+1
615 ilev2 = ilev1 - 1 + column%varNumLev(vnl_varListIndex(varName_opt))
616 onecolumn => column%all(ilev1:ilev2,headerIndex)
617 else
618 call utl_abort('col_getColumn: Unknown variable name! ' // varName_opt)
619 end if
620 else
621 onecolumn => column%all(:,headerIndex)
622 end if
623
624 end function col_getColumn
625
626 !--------------------------------------------------------------------------
627 ! col_getElem
628 !--------------------------------------------------------------------------
629 function col_getElem(column,ilev,headerIndex,varName_opt) result(value)
630 implicit none
631
632 ! Arguments:
633 type(struct_columnData), intent(in) :: column
634 integer, intent(in) :: ilev
635 integer, intent(in) :: headerIndex
636 character(len=*), optional, intent(in) :: varName_opt
637 ! Result:
638 real(8) :: value
639
640 if(present(varName_opt)) then
641 if(.not.col_varExist(column,varName_opt)) call utl_Abort('col_getElem: Unknown variable name! ' // varName_opt)
642 value = column%all(column%varOffset(vnl_varListIndex(varName_opt))+ilev,headerIndex)
643 else
644 value = column%all(ilev,headerIndex)
645 end if
646
647 end function col_getElem
648
649 !--------------------------------------------------------------------------
650 ! col_getNumLev
651 !--------------------------------------------------------------------------
652 function col_getNumLev(column,varLevel,varName_opt) result(nlev)
653 implicit none
654
655 ! Arguments:
656 type(struct_columnData), intent(in) :: column
657 character(len=*), intent(in) :: varLevel
658 character(len=*), optional, intent(in) :: varName_opt
659 ! Result:
660 integer :: nlev
661
662 nlev = vco_getNumLev(column%vco,varLevel,varName_opt)
663
664 end function col_getNumLev
665
666 !--------------------------------------------------------------------------
667 ! col_getNumCol
668 !--------------------------------------------------------------------------
669 function col_getNumCol(column) result(numColumn)
670 implicit none
671
672 ! Arguments:
673 type(struct_columnData), intent(in) :: column
674 ! Result:
675 integer :: numColumn
676
677 numColumn = column%numCol
678
679 end function col_getNumCol
680
681 !--------------------------------------------------------------------------
682 ! col_getVco
683 !--------------------------------------------------------------------------
684 function col_getVco(column) result(vco_ptr)
685 implicit none
686
687 ! Arguments:
688 type(struct_columnData), intent(in) :: column
689 ! Result:
690 type(struct_vco), pointer :: vco_ptr
691
692 vco_ptr => column%vco
693
694 end function col_getVco
695
696 !--------------------------------------------------------------------------
697 ! col_setVco
698 !--------------------------------------------------------------------------
699 subroutine col_setVco(column,vco_ptr)
700 implicit none
701
702 ! Arguments:
703 type(struct_columnData), intent(inout) :: column
704 type(struct_vco), pointer, intent(in) :: vco_ptr
705
706 column%vco => vco_ptr
707
708 end subroutine col_setVco
709
710end module columnData_mod