File: mytest.f90

package info (click to toggle)
cmor 3.14.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 16,976 kB
  • sloc: ansic: 28,053; f90: 13,893; python: 12,699; sh: 3,739; makefile: 111
file content (165 lines) | stat: -rw-r--r-- 4,843 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
program main

  USE cmor_users_functions
  implicit none

  integer ncid

  type dims
     integer n
     character(256) name
     character(256) units
     double precision, DIMENSION(:), pointer :: values
     double precision, DIMENSION(:,:), pointer :: bounds     
     type(dims), pointer :: next
  end type dims
  character(256) filein
  type(dims), pointer :: mydims,current
  integer ndim,i,j,ntot
!  integer, allocatable, dimension(:):: arrayin
  double precision, allocatable, dimension(:):: arrayin
  integer, dimension(7):: dimlength = (/ (1,i=1,7) /)
  integer, PARAMETER::verbosity = 2
  integer ierr
  integer, allocatable, dimension(:) :: myaxis
  integer myvar

  print*, 'hi enter test file case'
  filein='Test/tas.asc'
  read(5,'(A)') filein
  open(unit=23,file=filein,form='formatted') 
  call allocate_dims(23,mydims,ndim)
  allocate(myaxis(ndim))
  current=>mydims
  ntot=1
  do i =1,ndim
     ntot=ntot*current%n
     current=>current%next
  enddo
  call read_ascii(23,mydims, ndim,ntot,arrayin)
  
  
!!$!! Ok here is the part where we define or variable/axis,etc... 
!!$!! Assuming that Karl's code is ok...
!!$
  
  print*,'CMOR SETUP'
!!$  
  ierr = cmor_setup(inpath='Test',   &
       netcdf_file_action='replace',                                       &
       set_verbosity=1,                                                    &
       exit_control=1)
    
  print*,'CMOR DATASET'
  ierr = cmor_dataset_json("Test/CMOR_input_example.json")
  
  current=>mydims
  do i = 0,ndim-1
     print*,'CMOR AXIS',i,'AAAAAAA**************************'
     print*, 'Name:',trim(adjustl(current%name)),'--',trim(adjustl(mydims%name))
     print*, current%units
     print*, current%n,size(current%values)
     print*, current%values(1:min(4,size(current%values)))
     print*, current%bounds(1:2,1:min(4,size(current%values)))
     if (i==0) then
        myaxis(ndim-i)=cmor_axis('CMOR_SAMPLE_TABLE', &
             table_entry=current%name,&
             units=current%units,&
             length=current%n)
     else
        myaxis(ndim-i)=cmor_axis('CMOR_SAMPLE_TABLE', &
             table_entry=current%name,&
             units=current%units,&
             length=current%n,&
             coord_vals=current%values,&
             cell_bounds=current%bounds)
     endif
     current=>current%next
  enddo

  print*,'CMOR VAR'
  myvar=cmor_variable('CMOR_SAMPLE_TABLE',&
       'tas',&
       'K',&
       myaxis,&
       missing_value=1.e20)

!! figures out length of dimension other than time



  j=ntot/mydims%n
  print*, '&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&'
  do i=1,mydims%n !! write times one at a time
!!$     print*,'before:', arrayin(j*(i-1)+1),j*(i-1)+1,j,i,ntot
     print*, 'size of arrayin',size(arrayin(j*(i-1)+1:j*i))
     ierr = cmor_write( &
          var_id        = myvar,                        &
          data          = arrayin(j*(i-1)+1:j*i), &
          ntimes_passed = 1,                              &
          time_vals     = mydims%values(i:i),                         &
          time_bnds     = mydims%bounds(1:2,i:i) &
          )
  enddo

  ierr = cmor_close()  

 
contains
  subroutine allocate_dims(file_id,mydims,ndim)
    implicit none
    integer i,n,j,tmp,file_id
    integer, intent(inout)::ndim
    type(dims) , pointer :: tmpdims,mydims
    read(file_id,'(I4)') ndim
!!$    allocate(dimlength(ndim))
    n=1
    allocate(mydims)
    tmpdims=>mydims
    do i = 1, ndim
       read(file_id,'(I4)') tmp
       dimlength(8-i)=tmp
       allocate(tmpdims%values(tmp))
       allocate(tmpdims%bounds(2,tmp))
       tmpdims%n=tmp
       allocate(tmpdims%next)
       tmpdims=>tmpdims%next
       n=n*tmp
    enddo
    deallocate(tmpdims)
    allocate(arrayin(n))
  end subroutine allocate_dims
  
  subroutine read_ascii(file_unit,mydims,ndim,ntot,arrayin)
    implicit none
    type(dims), pointer::  mydims
!    integer, dimension(ntot),intent(inout) :: arrayin
    double precision, dimension(ntot),intent(inout) :: arrayin
    type(dims), pointer ::  current
    integer, intent(in)::ndim,file_unit
    integer n,ntot,i,j,k
    
    current=>mydims
    ntot=1
    do i =1,ndim
       n=current%n
       ntot=ntot*n
       read(file_unit,*) current%name
       print*, 'NAME is:',current%name,trim(adjustl(mydims%name)),n,ntot
       read(file_unit,'(A)') current%units
       read(file_unit,*) (current%values(j),j=1,n)
       read(file_unit,*) ((current%bounds(j,k),j=1,2),k=1,n)
       print*, current%values(1),current%values(n)
       print*, current%bounds(1,1),current%bounds(2,1)
       print*, current%bounds(1,n),current%bounds(2,n)
       current=>current%next
    enddo
    print *, 'ntot:',ntot
       read(file_unit,*) (arrayin(i),i=1,ntot)

print* ,trim(adjustl(mydims%name))
print*,'done reading'
  end subroutine read_ascii

end program main