source: branches/flexpart91_hasod/src_parallel/readcommand.f90 @ 10

Last change on this file since 10 was 10, checked in by hasod, 11 years ago

ADD: namelist input implemented for all common input files

File size: 21.8 KB
RevLine 
[8]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 readcommand
23
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine reads the user specifications for the current model run.  *
27  !                                                                            *
28  !     Author: A. Stohl                                                       *
29  !     18 May 1996                                                            *
[10]30  !     HSO, 14 August 2013
31  !     Added optional namelist input
[8]32  !                                                                            *
33  !*****************************************************************************
34  !                                                                            *
35  ! Variables:                                                                 *
36  ! bdate                beginning date as Julian date                         *
37  ! ctl                  factor by which time step must be smaller than        *
38  !                      Lagrangian time scale                                 *
39  ! ibdate,ibtime        beginnning date and time (YYYYMMDD, HHMISS)           *
40  ! ideltas [s]          modelling period                                      *
41  ! iedate,ietime        ending date and time (YYYYMMDD, HHMISS)               *
42  ! ifine                reduction factor for vertical wind time step          *
43  ! outputforeachrel     for forward runs it is possible either to create      *
44  !                      one outputfield or several for each releasepoint      *
45  ! iflux                switch to turn on (1)/off (0) flux calculations       *
46  ! iout                 1 for conc. (residence time for backward runs) output,*
47  !                      2 for mixing ratio output, 3 both, 4 for plume        *
48  !                      trajectory output, 5 = options 1 and 4                *
49  ! ipin                 1 continue simulation with dumped particle data, 0 no *
50  ! ipout                0 no particle dump, 1 every output time, 3 only at end*
51  ! itsplit [s]          time constant for particle splitting                  *
52  ! loutaver [s]         concentration output is an average over loutaver      *
53  !                      seconds                                               *
54  ! loutsample [s]       average is computed from samples taken every [s]      *
55  !                      seconds                                               *
56  ! loutstep [s]         time interval of concentration output                 *
57  ! lsynctime [s]        synchronisation time interval for all particles       *
58  ! lagespectra          switch to turn on (1)/off (0) calculation of age      *
59  !                      spectra                                               *
60  ! lconvection          value of either 0 and 1 indicating mixing by          *
61  !                      convection                                            *
62  !                      = 0 .. no convection                                  *
63  !                      + 1 .. parameterisation of mixing by subgrid-scale    *
64  !                              convection = on                               *
65  ! lsubgrid             switch to turn on (1)/off (0) subgrid topography      *
66  !                      parameterization                                      *
67  ! method               method used to compute the particle pseudovelocities  *
68  ! mdomainfill          1 use domain-filling option, 0 not, 2 use strat. O3   *
69  !                                                                            *
70  ! Constants:                                                                 *
71  ! unitcommand          unit connected to file COMMAND                        *
72  !                                                                            *
73  !*****************************************************************************
74
75  use par_mod
76  use com_mod
77
78  implicit none
79
80  real(kind=dp) :: juldate
81  character(len=50) :: line
82  logical :: old
83  logical :: nmlout=.false.
84  integer :: readerror
85
[10]86  ! declaration of namelist
[8]87  namelist /command/ &
88    ldirect, &
89    ibdate,ibtime, &
90    iedate,ietime, &
91    loutstep, &
92    loutaver, &
93    loutsample, &
94    itsplit, &
95    lsynctime, &
96    ctl, &
97    ifine, &
98    iout, &
99    ipout, &
[10]100    lnetcdfout, &
[8]101    lsubgrid, &
102    lconvection, &
103    lagespectra, &
104    ipin, &
105    ioutputforeachrelease, &
106    iflux, &
107    mdomainfill, &
108    ind_source, &
109    ind_receptor, &
110    mquasilag, &
111    nested_output, &
112    linit_cond
113
[10]114  ! Presetting the namelist command
[8]115  ldirect=0
[10]116  ibdate=0
[8]117  ibtime=0
[10]118  iedate=0
[8]119  ietime=0
120  loutstep=10800
121  loutaver=10800
122  loutsample=900
123  itsplit=999999999
124  lsynctime=900
125  ctl=-5.0
126  ifine=4
127  iout=3
128  ipout=0
[10]129  lnetcdfout=0
[8]130  lsubgrid=1
131  lconvection=1
132  lagespectra=0
133  ipin=1
134  ioutputforeachrelease=0
135  iflux=1
136  mdomainfill=0
137  ind_source=1
138  ind_receptor=1
139  mquasilag=0
140  nested_output=0
141  linit_cond=0
142
143  ! Open the command file and read user options
144  ! Namelist input first: try to read as namelist file
145  !**************************************************************************
146  open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old', &
[10]147         err=999)
[8]148
[10]149  ! try namelist input
[8]150  read(unitcommand,command,iostat=readerror)
151
152  ! If error in namelist format, try to open with old input code
153  ! Failsafe detection with ldirect initial value change
[10]154  if ((readerror.ne.0).or.(ldirect.eq.0)) then
[8]155
[10]156    rewind(unitcommand)
[8]157
158    ! Check the format of the COMMAND file (either in free format,
159    ! or using formatted mask)
160    ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
161    !**************************************************************************
162
163    call skplin(9,unitcommand)
164    read (unitcommand,901) line
165  901   format (a)
166    if (index(line,'LDIRECT') .eq. 0) then
167      old = .false.
168    else
169      old = .true.
170    endif
171    rewind(unitcommand)
172
173    ! Read parameters
174    !****************
175
176    call skplin(7,unitcommand)
177    if (old) call skplin(1,unitcommand)
178
179    read(unitcommand,*) ldirect
180    if (old) call skplin(3,unitcommand)
181    read(unitcommand,*) ibdate,ibtime
182    if (old) call skplin(3,unitcommand)
183    read(unitcommand,*) iedate,ietime
184    if (old) call skplin(3,unitcommand)
185    read(unitcommand,*) loutstep
186    if (old) call skplin(3,unitcommand)
187    read(unitcommand,*) loutaver
188    if (old) call skplin(3,unitcommand)
189    read(unitcommand,*) loutsample
190    if (old) call skplin(3,unitcommand)
191    read(unitcommand,*) itsplit
192    if (old) call skplin(3,unitcommand)
193    read(unitcommand,*) lsynctime
194    if (old) call skplin(3,unitcommand)
195    read(unitcommand,*) ctl
196    if (old) call skplin(3,unitcommand)
197    read(unitcommand,*) ifine
198    if (old) call skplin(3,unitcommand)
199    read(unitcommand,*) iout
200    if (old) call skplin(3,unitcommand)
201    read(unitcommand,*) ipout
202    if (old) call skplin(3,unitcommand)
203    read(unitcommand,*) lsubgrid
204    if (old) call skplin(3,unitcommand)
205    read(unitcommand,*) lconvection
206    if (old) call skplin(3,unitcommand)
207    read(unitcommand,*) lagespectra
208    if (old) call skplin(3,unitcommand)
209    read(unitcommand,*) ipin
210    if (old) call skplin(3,unitcommand)
211    read(unitcommand,*) ioutputforeachrelease
212    if (old) call skplin(3,unitcommand)
213    read(unitcommand,*) iflux
214    if (old) call skplin(3,unitcommand)
215    read(unitcommand,*) mdomainfill
216    if (old) call skplin(3,unitcommand)
217    read(unitcommand,*) ind_source
218    if (old) call skplin(3,unitcommand)
219    read(unitcommand,*) ind_receptor
220    if (old) call skplin(3,unitcommand)
221    read(unitcommand,*) mquasilag
222    if (old) call skplin(3,unitcommand)
223    read(unitcommand,*) nested_output
224    if (old) call skplin(3,unitcommand)
225    read(unitcommand,*) linit_cond
226    close(unitcommand)
227
228  endif ! input format
229
230  ! write command file in namelist format to output directory if requested
231  if (nmlout.eqv..true.) then
232    open(unitcommand,file=path(2)(1:length(2))//'COMMAND.namelist',status='new',err=999)
233    write(unitcommand,nml=command)
234    close(unitcommand)
235  endif
236
237  ifine=max(ifine,1)
238
239  ! Determine how Markov chain is formulated (for w or for w/sigw)
240  !***************************************************************
241
242  if (ctl.ge.0.1) then
243    turbswitch=.true.
244  else
245    turbswitch=.false.
246    ifine=1
247  endif
248  fine=1./real(ifine)
249  ctl=1./ctl
250
251  ! Set the switches required for the various options for input/output units
252  !*************************************************************************
253  !AF Set the switches IND_REL and IND_SAMP for the release and sampling
254  !Af switches for the releasefile:
255  !Af IND_REL =  1 : xmass * rho
256  !Af IND_REL =  0 : xmass * 1
257
258  !Af switches for the conccalcfile:
259  !AF IND_SAMP =  0 : xmass * 1
260  !Af IND_SAMP = -1 : xmass / rho
261
262  !AF IND_SOURCE switches between different units for concentrations at the source
263  !Af   NOTE that in backward simulations the release of computational particles
264  !Af   takes place at the "receptor" and the sampling of p[articles at the "source".
265  !Af          1 = mass units
266  !Af          2 = mass mixing ratio units
267  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
268  !Af          1 = mass units
269  !Af          2 = mass mixing ratio units
270
271  if ( ldirect .eq. 1 ) then  ! FWD-Run
272  !Af set release-switch
273     if (ind_source .eq. 1 ) then !mass
274        ind_rel = 0
275     else ! mass mix
276        ind_rel = 1
277     endif
278  !Af set sampling switch
279     if (ind_receptor .eq. 1) then !mass
280        ind_samp = 0
281     else ! mass mix
282        ind_samp = -1
283     endif
284  elseif (ldirect .eq. -1 ) then !BWD-Run
285  !Af set sampling switch
286     if (ind_source .eq. 1 ) then !mass
287        ind_samp = -1
288     else ! mass mix
289        ind_samp = 0
290     endif
291  !Af set release-switch
292     if (ind_receptor .eq. 1) then !mass
293        ind_rel = 1
294     else ! mass mix
295        ind_rel = 0
296     endif
297  endif
298
299  !*************************************************************
300  ! Check whether valid options have been chosen in file COMMAND
301  !*************************************************************
302
303  ! Check options for initial condition output: Switch off for forward runs
304  !************************************************************************
305
306  if (ldirect.eq.1) linit_cond=0
307  if ((linit_cond.lt.0).or.(linit_cond.gt.2)) then
308    write(*,*) ' #### FLEXPART MODEL ERROR! INVALID OPTION    #### '
309    write(*,*) ' #### FOR LINIT_COND IN FILE "COMMAND".       #### '
310    stop
311  endif
312
313  ! Check input dates
314  !******************
315
316  if (iedate.lt.ibdate) then
317    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING DATE    #### '
318    write(*,*) ' #### IS LARGER THAN ENDING DATE. CHANGE      #### '
319    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
320    write(*,*) ' #### "COMMAND".                              #### '
321    stop
322  else if (iedate.eq.ibdate) then
323    if (ietime.lt.ibtime) then
324    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING TIME    #### '
325    write(*,*) ' #### IS LARGER THAN ENDING TIME. CHANGE      #### '
326    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
327    write(*,*) ' #### "COMMAND".                              #### '
328    stop
329    endif
330  endif
331
332
333  ! Determine kind of dispersion method
334  !************************************
335
336  if (ctl.gt.0.) then
337    method=1
338    mintime=minstep
339  else
340    method=0
341    mintime=lsynctime
342  endif
343
[10]344!  check for netcdf output switch (use for non-namelist input only!)
[8]345  if (iout.ge.8) then
[10]346     lnetcdfout = 1
[8]347     iout = iout - 8
348#ifndef NETCDF_OUTPUT
[10]349     print*,'ERROR: netcdf output not activated during compile time but used in COMMAND file!'
[8]350#endif
351  endif
352
353  ! Check whether a valid option for gridded model output has been chosen
354  !**********************************************************************
355
356  if ((iout.lt.1).or.(iout.gt.5)) then
357    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
358    write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4, OR 5 FOR       #### '
359    write(*,*) ' #### STANDARD FLEXPART OUTPUT OR  9 - 13     #### '
360    write(*,*) ' #### FOR NETCDF OUTPUT                       #### '
361    stop
362  endif
363
364  !AF check consistency between units and volume mixing ratio
365  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
366       (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
367    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
368    write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED      #### '
369    write(*,*) ' #### FOR MASS UNITS (at the moment)          #### '
370    stop
371  endif
372
373
374
375  ! For quasilag output for each release is forbidden
376  !*****************************************************************************
377
378  if ((ioutputforeachrelease.eq.1).and.(mquasilag.eq.1)) then
379      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
380      write(*,*) '#### OUTPUTFOREACHRELEASE AND QUASILAGRANGIAN####'
381      write(*,*) '#### MODE IS NOT POSSIBLE   !                ####'
382      stop
383  endif
384
385
386  ! For quasilag backward is forbidden
387  !*****************************************************************************
388
389  if ((ldirect.lt.0).and.(mquasilag.eq.1)) then
390      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
391      write(*,*) '#### FOR BACKWARD RUNS, QUASILAGRANGIAN MODE ####'
392      write(*,*) '#### IS NOT POSSIBLE   !                     ####'
393      stop
394  endif
395
396
397  ! For backward runs one releasefield for all releases makes no sense,
398  ! For quasilag and domainfill ioutputforechrelease is forbidden
399  !*****************************************************************************
400
401  if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
402      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
403      write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
404      write(*,*) '#### MUST BE SET TO ONE!                     ####'
405      stop
406  endif
407
408
409  ! For backward runs one releasefield for all releases makes no sense,
410  ! and is "forbidden"
411  !*****************************************************************************
412
413  if ((mdomainfill.eq.1).and.(ioutputforeachrelease.eq.1)) then
414      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
415      write(*,*) '#### FOR DOMAIN FILLING RUNS OUTPUT FOR      ####'
416      write(*,*) '#### EACH RELEASE IS FORBIDDEN !             ####'
417      stop
418  endif
419
420
421  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
422  ! For backward runs, only residence time output (iout=1) or plume trajectories (iout=4),
423  ! or both (iout=5) makes sense; other output options are "forbidden"
424  !*****************************************************************************
425
426  if (ldirect.lt.0) then
427    if ((iout.eq.2).or.(iout.eq.3)) then
428      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
429      write(*,*) '#### FOR BACKWARD RUNS, IOUT MUST BE 1,4,OR 5####'
430      stop
431    endif
432  endif
433
434
435  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
436  ! and is "forbidden"
437  !*****************************************************************************
438
439  if (mdomainfill.ge.1) then
440    if ((iout.eq.4).or.(iout.eq.5)) then
441      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
442      write(*,*) '#### FOR DOMAIN-FILLING TRAJECTORY OPTION,   ####'
443      write(*,*) '#### IOUT MUST NOT BE SET TO 4 OR 5.         ####'
444      stop
445    endif
446  endif
447
448
449
450  ! Check whether a valid options for particle dump has been chosen
451  !****************************************************************
452
453  if ((ipout.ne.0).and.(ipout.ne.1).and.(ipout.ne.2)) then
454    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
455    write(*,*) ' #### IPOUT MUST BE 1, 2 OR 3!                #### '
456    stop
457  endif
458
459  if(lsubgrid.ne.1) then
460    write(*,*) '             ----------------               '
461    write(*,*) ' INFORMATION: SUBGRIDSCALE TERRAIN EFFECT IS'
462    write(*,*) ' NOT PARAMETERIZED DURING THIS SIMULATION.  '
463    write(*,*) '             ----------------               '
464  endif
465
466
467  ! Check whether convection scheme is either turned on or off
468  !***********************************************************
469
470  if ((lconvection.ne.0).and.(lconvection.ne.1)) then
471    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
472    write(*,*) ' #### LCONVECTION MUST BE SET TO EITHER 1 OR 0#### '
473    stop
474  endif
475
476
477  ! Check whether synchronisation interval is sufficiently short
478  !*************************************************************
479
480  if (lsynctime.gt.(idiffnorm/2)) then
481    write(*,*) ' #### FLEXPART MODEL ERROR! SYNCHRONISATION   #### '
482    write(*,*) ' #### TIME IS TOO LONG. MAKE IT SHORTER.      #### '
483    write(*,*) ' #### MINIMUM HAS TO BE: ', idiffnorm/2
484    stop
485  endif
486
487
488  ! Check consistency of the intervals, sampling periods, etc., for model output
489  !*****************************************************************************
490
491  if (loutaver.eq.0) then
492    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
493    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
494    write(*,*) ' #### ZERO.                                   #### '
495    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
496    stop
497  endif
498
499  if (loutaver.gt.loutstep) then
500    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
501    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
502    write(*,*) ' #### GREATER THAN INTERVAL OF OUTPUT.        #### '
503    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
504    stop
505  endif
506
507  if (loutsample.gt.loutaver) then
508    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
509    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
510    write(*,*) ' #### GREATER THAN TIME AVERAGE OF OUTPUT.    #### '
511    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
512    stop
513  endif
514
515  if (mod(loutaver,lsynctime).ne.0) then
516    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
517    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
518    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
519    stop
520  endif
521
522  if ((loutaver/lsynctime).lt.2) then
523    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
524    write(*,*) ' #### CONCENTRATION FIELD MUST BE AT LEAST    #### '
525    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
526    stop
527  endif
528
529  if (mod(loutstep,lsynctime).ne.0) then
530    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
531    write(*,*) ' #### CONCENTRATION FIELDS MUST BE A MULTIPLE #### '
532    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
533    stop
534  endif
535
536  if ((loutstep/lsynctime).lt.2) then
537    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
538    write(*,*) ' #### CONCENTRATION FIELDS MUST BE AT LEAST   #### '
539    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
540    stop
541  endif
542
543  if (mod(loutsample,lsynctime).ne.0) then
544    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
545    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
546    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
547    stop
548  endif
549
550  if (itsplit.lt.loutaver) then
551    write(*,*) ' #### FLEXPART MODEL ERROR! SPLITTING TIME FOR#### '
552    write(*,*) ' #### PARTICLES IS TOO SHORT. PLEASE INCREASE #### '
553    write(*,*) ' #### SPLITTING TIME CONSTANT.                #### '
554    stop
555  endif
556
557  if ((mquasilag.eq.1).and.(iout.ge.4)) then
558    write(*,*) ' #### FLEXPART MODEL ERROR! CONFLICTING       #### '
559    write(*,*) ' #### OPTIONS: IF MQUASILAG=1, PLUME          #### '
560    write(*,*) ' #### TRAJECTORY OUTPUT IS IMPOSSIBLE.        #### '
561    stop
562  endif
563
564  ! Compute modeling time in seconds and beginning date in Julian date
565  !*******************************************************************
566
567  outstep=real(abs(loutstep))
568  if (ldirect.eq.1) then
569    bdate=juldate(ibdate,ibtime)
570    edate=juldate(iedate,ietime)
571    ideltas=nint((edate-bdate)*86400.)
572  else if (ldirect.eq.-1) then
573    loutaver=-1*loutaver
574    loutstep=-1*loutstep
575    loutsample=-1*loutsample
576    lsynctime=-1*lsynctime
577    bdate=juldate(iedate,ietime)
578    edate=juldate(ibdate,ibtime)
579    ideltas=nint((edate-bdate)*86400.)
580  else
581    write(*,*) ' #### FLEXPART MODEL ERROR! DIRECTION IN      #### '
582    write(*,*) ' #### FILE "COMMAND" MUST BE EITHER -1 OR 1.  #### '
583    stop
584  endif
585
586  return
587
588999   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND"    #### '
589  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
590  write(*,'(a)') path(1)(1:length(1))
591  stop
592
593end subroutine readcommand
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG