source: branches/sabine/readwind.f90 @ 6

Last change on this file since 6 was 6, checked in by saeck, 11 years ago

import to sabine

File size: 18.8 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine readwind(indj,n,uuh,vvh,wwh)
23
24  !**********************************************************************
25  !                                                                     *
26  !             TRAJECTORY MODEL SUBROUTINE READWIND                    *
27  !                                                                     *
28  !**********************************************************************
29  !                                                                     *
30  !             AUTHOR:      G. WOTAWA                                  *
31  !             DATE:        1997-08-05                                 *
32  !             LAST UPDATE: 2000-10-17, Andreas Stohl                  *
33  !             CHANGE: 11/01/2008, Harald Sodemann, GRIB1/2 input with *
34  !                                 ECMWF grib_api                      *
35  !             CHANGE: 03/12/2008, Harald Sodemann, update to f90 with *
36  !                                 ECMWF grib_api                      *
37  !                                                                     *
38  !**********************************************************************
39  !  Changes, Bernd C. Krueger, Feb. 2001:
40  !   Variables tth and qvh (on eta coordinates) in common block
41  !**********************************************************************
42  !                                                                     *
43  ! DESCRIPTION:                                                        *
44  !                                                                     *
45  ! READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE   *
46  ! INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE          *
47  !                                                                     *
48  ! INPUT:                                                              *
49  ! indj               indicates number of the wind field to be read in *
50  ! n                  temporal index for meteorological fields (1 to 3)*
51  !                                                                     *
52  ! IMPORTANT VARIABLES FROM COMMON BLOCK:                              *
53  !                                                                     *
54  ! wfname             File name of data to be read in                  *
55  ! nx,ny,nuvz,nwz     expected field dimensions                        *
56  ! nlev_ec            number of vertical levels ecmwf model            *
57  ! uu,vv,ww           wind fields                                      *
58  ! tt,qv              temperature and specific humidity                *
59  ! ps                 surface pressure                                 *
60  !                                                                     *
61  !**********************************************************************
62
63  use GRIB_API
64  use par_mod
65  use com_mod
66
67  implicit none
68
69  !HSO  parameters for grib_api
70  integer :: ifile
71  integer :: iret
72  integer :: igrib
73  integer :: gribVer,parCat,parNum,typSurf,valSurf,discipl
74  integer :: gotGrid
75  !HSO  end
76
77  real(kind=4) :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
78  real(kind=4) :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
79  real(kind=4) :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
80  integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax
81
82  ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING
83
84  ! dimension of isec2 at least (22+n), where n is the number of parallels or
85  ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid
86
87  ! dimension of zsec2 at least (10+nn), where nn is the number of vertical
88  ! coordinate parameters
89
90  integer :: isec1(56),isec2(22+nxmax+nymax)
91  real(kind=4) :: zsec4(jpunp)
92  real(kind=4) :: xaux,yaux,xaux0,yaux0
93  real(kind=8) :: xauxin,yauxin
94  real,parameter :: eps=1.e-4
95  real(kind=4) :: nsss(0:nxmax-1,0:nymax-1),ewss(0:nxmax-1,0:nymax-1)
96  real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1
97
98  logical :: hflswitch,strswitch
99
100  !HSO  grib api error messages
101  character(len=24) :: gribErrorMsg = 'Error reading grib file'
102  character(len=20) :: gribFunction = 'readwind'
103
104  hflswitch=.false.
105  strswitch=.false.
106  levdiff2=nlev_ec-nwz+1
107  iumax=0
108  iwmax=0
109
110  !
111  ! OPENING OF DATA FILE (GRIB CODE)
112  !
1135   call grib_open_file(ifile,path(3)(1:length(3)) &
114         //trim(wfname(indj)),'r',iret)
115  if (iret.ne.GRIB_SUCCESS) then
116    goto 888   ! ERROR DETECTED
117  endif
118  !turn on support for multi fields messages */
119  !call grib_multi_support_on
120
121  gotGrid=0
122  ifield=0
12310   ifield=ifield+1
124  !
125  ! GET NEXT FIELDS
126  !
127  call grib_new_from_file(ifile,igrib,iret)
128  if (iret.eq.GRIB_END_OF_FILE)  then
129    goto 50    ! EOF DETECTED
130  elseif (iret.ne.GRIB_SUCCESS) then
131    goto 888   ! ERROR DETECTED
132  endif
133
134  !first see if we read GRIB1 or GRIB2
135  call grib_get_int(igrib,'editionNumber',gribVer,iret)
136  call grib_check(iret,gribFunction,gribErrorMsg)
137
138  if (gribVer.eq.1) then ! GRIB Edition 1
139
140  !print*,'GRiB Edition 1'
141  !read the grib2 identifiers
142  call grib_get_int(igrib,'indicatorOfParameter',isec1(6),iret)
143  call grib_check(iret,gribFunction,gribErrorMsg)
144  call grib_get_int(igrib,'level',isec1(8),iret)
145  call grib_check(iret,gribFunction,gribErrorMsg)
146
147  !change code for etadot to code for omega
148  if (isec1(6).eq.77) then
149    isec1(6)=135
150  endif
151
152  else
153
154  !print*,'GRiB Edition 2'
155  !read the grib2 identifiers
156  call grib_get_int(igrib,'discipline',discipl,iret)
157  call grib_check(iret,gribFunction,gribErrorMsg)
158  call grib_get_int(igrib,'parameterCategory',parCat,iret)
159  call grib_check(iret,gribFunction,gribErrorMsg)
160  call grib_get_int(igrib,'parameterNumber',parNum,iret)
161  call grib_check(iret,gribFunction,gribErrorMsg)
162  call grib_get_int(igrib,'typeOfFirstFixedSurface',typSurf,iret)
163  call grib_check(iret,gribFunction,gribErrorMsg)
164  call grib_get_int(igrib,'level',valSurf,iret)
165  call grib_check(iret,gribFunction,gribErrorMsg)
166
167  !print*,discipl,parCat,parNum,typSurf,valSurf
168
169  !convert to grib1 identifiers
170  isec1(6)=-1
171  isec1(7)=-1
172  isec1(8)=-1
173  isec1(8)=valSurf     ! level
174  if ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! T
175    isec1(6)=130         ! indicatorOfParameter
176  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.105)) then ! U
177    isec1(6)=131         ! indicatorOfParameter
178  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.105)) then ! V
179    isec1(6)=132         ! indicatorOfParameter
180  elseif ((parCat.eq.1).and.(parNum.eq.0).and.(typSurf.eq.105)) then ! Q
181    isec1(6)=133         ! indicatorOfParameter
182  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.1)) then !SP
183    isec1(6)=134         ! indicatorOfParameter
184  elseif ((parCat.eq.2).and.(parNum.eq.32)) then ! W, actually eta dot
185    isec1(6)=135         ! indicatorOfParameter
186  elseif ((parCat.eq.3).and.(parNum.eq.0).and.(typSurf.eq.101)) then !SLP
187    isec1(6)=151         ! indicatorOfParameter
188  elseif ((parCat.eq.2).and.(parNum.eq.2).and.(typSurf.eq.103)) then ! 10U
189    isec1(6)=165         ! indicatorOfParameter
190  elseif ((parCat.eq.2).and.(parNum.eq.3).and.(typSurf.eq.103)) then ! 10V
191    isec1(6)=166         ! indicatorOfParameter
192  elseif ((parCat.eq.0).and.(parNum.eq.0).and.(typSurf.eq.103)) then ! 2T
193    isec1(6)=167         ! indicatorOfParameter
194  elseif ((parCat.eq.0).and.(parNum.eq.6).and.(typSurf.eq.103)) then ! 2D
195    isec1(6)=168         ! indicatorOfParameter
196  elseif ((parCat.eq.1).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SD
197    isec1(6)=141         ! indicatorOfParameter
198  elseif ((parCat.eq.6).and.(parNum.eq.1)) then ! CC
199    isec1(6)=164         ! indicatorOfParameter
200  elseif ((parCat.eq.1).and.(parNum.eq.9)) then ! LSP
201    isec1(6)=142         ! indicatorOfParameter
202  elseif ((parCat.eq.1).and.(parNum.eq.10)) then ! CP
203    isec1(6)=143         ! indicatorOfParameter
204  elseif ((parCat.eq.0).and.(parNum.eq.11).and.(typSurf.eq.1)) then ! SHF
205    isec1(6)=146         ! indicatorOfParameter
206  elseif ((parCat.eq.4).and.(parNum.eq.9).and.(typSurf.eq.1)) then ! SR
207    isec1(6)=176         ! indicatorOfParameter
208  elseif ((parCat.eq.2).and.(parNum.eq.17)) then ! EWSS
209    isec1(6)=180         ! indicatorOfParameter
210  elseif ((parCat.eq.2).and.(parNum.eq.18)) then ! NSSS
211    isec1(6)=181         ! indicatorOfParameter
212  elseif ((parCat.eq.3).and.(parNum.eq.4)) then ! ORO
213    isec1(6)=129         ! indicatorOfParameter
214  elseif ((parCat.eq.3).and.(parNum.eq.7)) then ! SDO
215    isec1(6)=160         ! indicatorOfParameter
216  elseif ((discipl.eq.2).and.(parCat.eq.0).and.(parNum.eq.0).and. &
217       (typSurf.eq.1)) then ! LSM
218    isec1(6)=172         ! indicatorOfParameter
219  else
220    print*,'***ERROR: undefined GRiB2 message found!',discipl, &
221         parCat,parNum,typSurf
222  endif
223
224  endif
225
226  !HSO  get the size and data of the values array
227  if (isec1(6).ne.-1) then
228    call grib_get_real4_array(igrib,'values',zsec4,iret)
229    call grib_check(iret,gribFunction,gribErrorMsg)
230  endif
231
232  !HSO  get the required fields from section 2 in a gribex compatible manner
233  if (ifield.eq.1) then
234  call grib_get_int(igrib,'numberOfPointsAlongAParallel', &
235       isec2(2),iret)
236  call grib_check(iret,gribFunction,gribErrorMsg)
237  call grib_get_int(igrib,'numberOfPointsAlongAMeridian', &
238       isec2(3),iret)
239  call grib_check(iret,gribFunction,gribErrorMsg)
240  call grib_get_int(igrib,'numberOfVerticalCoordinateValues', &
241       isec2(12))
242  call grib_check(iret,gribFunction,gribErrorMsg)
243  ! CHECK GRID SPECIFICATIONS
244  if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT'
245  if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT'
246  if(isec2(12)/2-1.ne.nlev_ec) &
247       stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT'
248  endif ! ifield
249
250  !HSO  get the second part of the grid dimensions only from GRiB1 messages
251  if ((gribVer.eq.1).and.(gotGrid.eq.0)) then
252    call grib_get_real8(igrib,'longitudeOfFirstGridPointInDegrees', &
253         xauxin,iret)
254    call grib_check(iret,gribFunction,gribErrorMsg)
255    call grib_get_real8(igrib,'latitudeOfLastGridPointInDegrees', &
256         yauxin,iret)
257    call grib_check(iret,gribFunction,gribErrorMsg)
258    xaux=xauxin+real(nxshift)*dx
259    yaux=yauxin
260    xaux0=xlon0
261    yaux0=ylat0
262    if(xaux.lt.0.) xaux=xaux+360.
263    if(yaux.lt.0.) yaux=yaux+360.
264    if(xaux0.lt.0.) xaux0=xaux0+360.
265    if(yaux0.lt.0.) yaux0=yaux0+360.
266    if(abs(xaux-xaux0).gt.eps) &
267         stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT'
268    if(abs(yaux-yaux0).gt.eps) &
269         stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT'
270    gotGrid=1
271  endif ! gotGrid
272
273  do j=0,nymin1
274    do i=0,nxfield-1
275      k=isec1(8)
276      if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE
277           zsec4(nxfield*(ny-j-1)+i+1)
278      if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY
279           zsec4(nxfield*(ny-j-1)+i+1)
280      if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY
281           zsec4(nxfield*(ny-j-1)+i+1)
282      if(isec1(6).eq.133) then                      !! SPEC. HUMIDITY
283        qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1)
284        if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) &
285             qvh(i,j,nlev_ec-k+2,n) = 0.
286  !        this is necessary because the gridded data may contain
287  !        spurious negative values
288      endif
289      if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS.
290           zsec4(nxfield*(ny-j-1)+i+1)
291
292      if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY
293           zsec4(nxfield*(ny-j-1)+i+1)
294      if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH
295           zsec4(nxfield*(ny-j-1)+i+1)
296      if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS.
297           zsec4(nxfield*(ny-j-1)+i+1)
298      if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER
299           zsec4(nxfield*(ny-j-1)+i+1)
300      if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY
301           zsec4(nxfield*(ny-j-1)+i+1)
302      if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY
303           zsec4(nxfield*(ny-j-1)+i+1)
304      if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE
305           zsec4(nxfield*(ny-j-1)+i+1)
306      if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT
307           zsec4(nxfield*(ny-j-1)+i+1)
308      if(isec1(6).eq.142) then                      !! LARGE SCALE PREC.
309        lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
310        if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0.
311      endif
312      if(isec1(6).eq.143) then                      !! CONVECTIVE PREC.
313        convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
314        if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0.
315      endif
316      if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX
317           zsec4(nxfield*(ny-j-1)+i+1)
318      if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) &
319           hflswitch=.true.    ! Heat flux available
320      if(isec1(6).eq.176) then                      !! SOLAR RADIATION
321        ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1)
322        if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0.
323      endif
324      if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS
325           zsec4(nxfield*(ny-j-1)+i+1)
326      if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS
327           zsec4(nxfield*(ny-j-1)+i+1)
328      if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. &
329           (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true.    ! stress available
330  !sec        strswitch=.true.
331      if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY
332           zsec4(nxfield*(ny-j-1)+i+1)/ga
333      if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY
334           zsec4(nxfield*(ny-j-1)+i+1)
335      if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK
336           zsec4(nxfield*(ny-j-1)+i+1)
337      if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1)
338      if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1)
339
340    end do
341  end do
342
343  call grib_release(igrib)
344  goto 10                      !! READ NEXT LEVEL OR PARAMETER
345  !
346  ! CLOSING OF INPUT DATA FILE
347  !
348
34950   call grib_close_file(ifile)
350
351  !error message if no fields found with correct first longitude in it
352  if (gotGrid.eq.0) then
353    print*,'***ERROR: input file needs to contain GRiB1 formatted'// &
354         'messages'
355    stop
356  endif
357
358  if(levdiff2.eq.0) then
359    iwmax=nlev_ec+1
360    do i=0,nxmin1
361      do j=0,nymin1
362        wwh(i,j,nlev_ec+1)=0.
363      end do
364    end do
365  endif
366
367  ! For global fields, assign the leftmost data column also to the rightmost
368  ! data column; if required, shift whole grid by nxshift grid points
369  !*************************************************************************
370
371  if (xglobal) then
372    call shift_field_0(ewss,nxfield,ny)
373    call shift_field_0(nsss,nxfield,ny)
374    call shift_field_0(oro,nxfield,ny)
375    call shift_field_0(excessoro,nxfield,ny)
376    call shift_field_0(lsm,nxfield,ny)
377    call shift_field(ps,nxfield,ny,1,1,2,n)
378    call shift_field(sd,nxfield,ny,1,1,2,n)
379    call shift_field(msl,nxfield,ny,1,1,2,n)
380    call shift_field(tcc,nxfield,ny,1,1,2,n)
381    call shift_field(u10,nxfield,ny,1,1,2,n)
382    call shift_field(v10,nxfield,ny,1,1,2,n)
383    call shift_field(tt2,nxfield,ny,1,1,2,n)
384    call shift_field(td2,nxfield,ny,1,1,2,n)
385    call shift_field(lsprec,nxfield,ny,1,1,2,n)
386    call shift_field(convprec,nxfield,ny,1,1,2,n)
387    call shift_field(sshf,nxfield,ny,1,1,2,n)
388    call shift_field(ssr,nxfield,ny,1,1,2,n)
389    call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n)
390    call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n)
391    call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1)
392    call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1)
393    call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1)
394  endif
395
396  do i=0,nxmin1
397    do j=0,nymin1
398      surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2)
399    end do
400  end do
401
402  if ((.not.hflswitch).or.(.not.strswitch)) then
403    write(*,*) 'WARNING: No flux data contained in GRIB file ', &
404         wfname(indj)
405
406  ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD
407  ! As ECMWF has increased the model resolution, such that now the first model
408  ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level
409  ! (3rd model level in FLEXPART) for the profile method
410  !***************************************************************************
411
412    do i=0,nxmin1
413      do j=0,nymin1
414        plev1=akz(3)+bkz(3)*ps(i,j,1,n)
415        pmean=0.5*(ps(i,j,1,n)+plev1)
416        tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n))
417        fu=-r_air*tv/ga/pmean
418        hlev1=fu*(plev1-ps(i,j,1,n))   ! HEIGTH OF FIRST MODEL LAYER
419        ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2)
420        fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2)
421        call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, &
422             tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, &
423             surfstr(i,j,1,n),sshf(i,j,1,n))
424        if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200.
425        if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400.
426      end do
427    end do
428  endif
429
430
431  ! Assign 10 m wind to model level at eta=1.0 to have one additional model
432  ! level at the ground
433  ! Specific humidity is taken the same as at one level above
434  ! Temperature is taken as 2 m temperature
435  !**************************************************************************
436
437     do i=0,nxmin1
438        do j=0,nymin1
439           uuh(i,j,1)=u10(i,j,1,n)
440           vvh(i,j,1)=v10(i,j,1,n)
441           qvh(i,j,1,n)=qvh(i,j,2,n)
442           tth(i,j,1,n)=tt2(i,j,1,n)
443        end do
444     end do
445
446  if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT'
447  if(iwmax.ne.nwz)    stop 'READWIND: NWZ NOT CONSISTENT'
448
449  return
450888   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
451  write(*,*) ' #### ',wfname(indj),'                    #### '
452  write(*,*) ' #### IS NOT GRIB FORMAT !!!                  #### '
453  stop 'Execution terminated'
454999   write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD         #### '
455  write(*,*) ' #### ',wfname(indj),'                    #### '
456  write(*,*) ' #### CANNOT BE OPENED !!!                    #### '
457  stop 'Execution terminated'
458
459end subroutine readwind
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG