#include "defs.h" ! module amie_module ! ! Module used to read data from the AMIE outputs (POT,mean energy, ! and energy flux). ! use params_module,only: nlon,nlat,nlev,nlonp1,nlonp4,nlatp1, | nlevp1,nmlon,nmlonp1,nmlat use addfld_module,only: addfld use input_module,only: iamie use aurora_module,only: gekev,gmwpm2 ! (nlonp1,0:nlatp1) in geog coords implicit none ! ! Define parameters for AMIE input data file: ! iamie in namelist: iamie=1 for Gang Lu HAO AMIE, iamie=2 for Geoff Crowley ASTRA AMIE, ! iamie=3 for Aaron Ridley U MI AMIE, iamie=4 for SWMF in U MI AMIE form ! colath.F restricts crit(1) between 15-30 colat so crit(2) is between 30 and 45) integer,parameter :: ! SWMF in 3 1-day files of 289 times from 0-0(24)UT ! | mxtimes = 289, ! maximum number of times of AMIE data per file (Dec 13,14,15 5 min) + 1 for 0UT ! | mxgdays = 1, ! maximum number of days of AMIE data (10 for Gang, 1 for Aaron and ASTRA) ! | lonmx = 180, ! maximum number of longitudes of AMIE data (1 less than grid), but SWMF is lonp1! ! | ithtrns = 89, ! corresponding to trans lat 1 mlat number of AMIE lats minus 1 ! | ithmx = 91, ! max number of AMIE mlats (90/(5/3)=54+1) to equator ! Geoff Crowley's amie files: ! ASTRA AMIE is 40min MLTS 00-24(dup) or 37 times, 31 lats 90-40 1.67 deg, every 5 min | mxtimes = 1153, ! maximum number of times of AMIE data per file (Dec 13-16 5 min) + 1 for 0UT Dec 17 | mxgdays = 1, ! maximum number of days of AMIE data (10 for Gang, 1 for Aaron and ASTRA) | lonmx = 36, ! maximum number of longitudes of AMIE data (1 less than grid), but ASTRA is lonp1! | ithtrns = 30, ! corresponding to trans lat 40 mlat number of AMIE lats minus 1 | ithmx = 55, ! max number of AMIE mlats (90/(5/3)=54+1) to equator ! Aaron Ridley's amie files: ! Ridley AMIE is 2deg 90-44 or 24 lats + 15 lats extended to 14 mlat where potential is zero ! | mxtimes = 5881, ! maximum number of times of AMIE data per file (Dec 13-16 1min) ! | mxgdays = 1, ! maximum number of days of AMIE data (10 for Gang, 1 for Aaron) ! | lonmx = 24, ! maximum number of longitudes of AMIE data (1 less than grid) ! | ithtrns = 24, ! Aaron: corresponding to trans lat 42 mlat (ie, 90 to 44 dmlat=2) ! | ithmx = 46, ! max number of AMIE mlats (90/2=45+1) to equator ! Gang Lu's AMIE files: ! | mxgdays = 10, ! maximum number of days of AMIE data ! | lonmx = 36, ! maximum number of longitudes of AMIE data ! | ithtrns = 30, ! corresponding to trans lat 40-deg ! | ithmx = 55, ! maximum number of latitudes of AMIE data ! Both set jmxm the same | jmxm = 2*ithmx-1 ! maximum number of global latitudes (only one 0mlat) ! lonp1 and latp1 are AMIE dim of lon (37 or lonmx+1) and lat (31 or ithmx+1 where 30 is 90 to 45.67 in dmlat=1.67) integer :: lonp1,latp1, | year(mxtimes),month(mxtimes),day(mxtimes),jday(mxtimes) ! Define AMIE output fields real :: ! | tiepot(nmlonp1,nmlat), | tieefx(nmlonp1,nmlat),tieekv(nmlonp1,nmlat) ! | ,efxg(nlonp1,0:nlatp1),ekvg(nlonp1,0:nlatp1) ! ! Define fields for AMIE input data file: ! electric potential in Volt ! mean energy in KeV ! energy flux in W/m^2 (wrong for U MI AMIE - is mW/m2 or ergs/cm2-s) ! amie_hpi_nh(sh) are AMIE hemi-integrated power real :: | cusplat_nh_amie, cuspmlt_nh_amie, cusplat_sh_amie, | cuspmlt_sh_amie, hpi_sh_amie, hpi_nh_amie, pcp_sh_amie, | pcp_nh_amie, crad(2), phida(2) ! ! 6/12: added by Emery for binary Ridley AMIE real :: Re,areatot real :: pii,dtr,htr real :: del,xmlt,dmlat,dlatm,dlonm,dmltm,rot real, dimension(ithtrns+1) :: mlatdeg,da real, dimension(lonmx+1) :: mltnh real, dimension(lonmx+1,ithtrns*2+2) :: mlts, lats, | TempPotential, AveEOut, EfluxOut ! 4/13: added by Emery for ascii ASTRA AMIE real,dimension(mxtimes,2*ithtrns+2,lonmx+1),private :: | epkv,wpm2,ekev integer,dimension(mxtimes,3),private :: amie_dhm contains !----------------------------------------------------------------------- subroutine init_amie ! ! Was called from tgcm.F, but now called in getamie for istep<=1 ! (this is not in init.F to avoid circular dependencies) ! use input_module,only: amienh,amiesh integer,external :: nextlu ! ! 6/12 added by Emery for Ridley binary AMIE files character (len=100), dimension(100) :: Lines character (len=90) :: char90 integer :: iError, i,j, lu_amie,nt,ihem,ntime,jt integer :: iy,imo,ida,ih,imin,isec real :: mlt,mlat,epot,eflx,qj,energy ! 6/12: Set AMIE NH or SH grid dimensions lonp1 and latp1 (to go 90 to 40 mlat in dmlat=2 for Aaron) latp1 = ithtrns + 1 lonp1 = lonmx + 1 pii = 4.*atan(1.) dtr = pii/180. htr = pii/12. ! Re=radius of the earth in m Re = 6371.e+3 C **** SET GRID SPACING DLATM, DLONG, DLONM C DMLAT=lat spacing in degrees of AMIE apex grid C 6/12: from amie.nc.out in ~ganglu/amie, have 31 lats w dmlat=1.67 from 90 to 40 mlat C and 37 mlts w dmltm=0.67h from 0 to 24MLT, where lonmx=36 (NOT 37!) C 18*3=54 so 90/54=1.67=dmlat (ithtrns=30 to 40 mlat, ithmx=55 to equator, jmxm=2*ithmx-1) ! dtr = pi/180. DMLAT = 180. / FLOAT(jmxm-1) ! AMIE grid spacing in deg (should be 2 deg Aaron, 1.67 deg ASTRA) DLATM = DMLAT*dtr DLONM = 2.*pii/FLOAT(lonmx) DMLTM = 24./FLOAT(lonmx) ! AMIE grid spacing in MLT (should be 1h Aaron, 0.67h ASTRA) areatot = 0. write(*,*) "dmlat dmlon=",dmlat,dlonm/dtr do j=1,latp1 ! In 'mirror', the min/max potentials are flipped from NH to SH ! Want lats from 90 to 40 mlat in both hem mlatdeg(j) = 90.0 - (j-1)*dmlat ! cos(90deg)=0 so pole is ignored with da(1)=0 da(j) = Re*cos(mlatdeg(j)*dtr)*dlonm*(Re*dmlat*dtr) areatot = areatot + da(j)*lonmx do i=1,lonp1 lats(i,j+latp1) = 90.0 - (j-1)*dmlat lats(i,j) = -90.0 + (j-1)*dmlat mltnh(i) = (i-1)*dmltm mlts(i,j) = (i-1)*dmltm mlts(i,j+latp1) = (i-1)*dmltm enddo enddo write (*,*) "areatot da(j)=",areatot,da if (len_trim(amienh) > 0) then write(6,"('Reading AMIENH file ',a)") trim(amienh) Lines(2) = trim(amienh) endif ! 04/13: Do different things for U MI and ASTRA AMIE ! ASTRA AMIE: if (iamie==2) then lu_amie = nextlu() write (6,"(1x,'lu_amie=',i3)") lu_amie open(lu_amie,file=trim(amienh),status='old') ! Skip header line read (lu_amie,"(a)") char90 ! Read every 5 min 4 days of data over 0000to2400 mlts from 90to40, -40to-90 mlat ! Store in order SH -90 to ~-40 mlat, NH +90 to ~+40mlat ntime = 0 do nt=1,mxtimes-1 ! ASTRA has both MLTS=0 and 24 or i=1 and lonp1 which are dups do i=1,lonp1 do ihem=1,2 do j=1,latp1 read(lu_amie,"(6i4,2f8.1,1x,4f12.4)") iy,imo,ida,ih,imin, | isec,mlt,mlat,epot,eflx,qj,energy if (ihem==1 .and. j==1 .and. i==1) then amie_dhm(nt,1) = ida + 334 ! for Dec06 case amie_dhm(nt,2) = ih amie_dhm(nt,3) = imin ntime = ntime + 1 write (6,"(1x,'ntime dayn hr min = ',4i4)") ntime, | (amie_dhm(nt,jt),jt=1,3) endif if (ihem==1) then epkv(nt,j+latp1,i) = epot wpm2(nt,j+latp1,i) = eflx ekev(nt,j+latp1,i) = energy else epkv(nt,latp1+1-j,i) = epot wpm2(nt,latp1+1-j,i) = eflx ekev(nt,latp1+1-j,i) = energy endif enddo enddo enddo ! do i=1,lonp1 ! wrap-around loc 1=lonp1 (mlt=0 and 24) not necessary since MLT=24 is also read in. enddo ! do nt=1,mxtimes-1 ! Extrapolate to 24 UT on Dec 17 amie_dhm(mxtimes,1) = amie_dhm(1,1)+4 amie_dhm(mxtimes,2) = 0 amie_dhm(mxtimes,3) = 0 do i=1,lonmx+1 do j=1,latp1*2 epkv(mxtimes,j,i) = epkv(mxtimes-1,j,i) | + (epkv(mxtimes-1,j,i)-epkv(mxtimes-2,j,i)) wpm2(mxtimes,j,i) = wpm2(mxtimes-1,j,i) | + (wpm2(mxtimes-1,j,i)-wpm2(mxtimes-2,j,i)) ekev(mxtimes,j,i) = ekev(mxtimes-1,j,i) | + (ekev(mxtimes-1,j,i)-ekev(mxtimes-2,j,i)) enddo enddo write (6,"(1x,'Reading AMIE times from 1 to mxtimes (dn,h,m) =', | i4,2x,3i4,1x,3i4,1x,3i4)") mxtimes,(amie_dhm(1,i),i=1,3), | (amie_dhm(mxtimes-1,j),j=1,3),(amie_dhm(mxtimes,jt),jt=1,3) endif ! if (iamie==2) ! U MI AMIE: if (iamie>=3) then if (len_trim(amiesh) > 0) then write(6,"('Reading AMIESH file ',a)") trim(amiesh) Lines(3) = trim(amiesh) else Lines(3) = "mirror" endif Lines(1) = "#AMIEFILES" Lines(4) = "" Lines(5) = "" Lines(6) = "" Lines(7) = "#DEBUG" Lines(8) = "2" Lines(9) = "0" Lines(10) = "" Lines(11) = "#END" write(*,*) "Calling IE_set_inputs" ! call IE_set_inputs(Lines) write(*,*) "=> Setting up Ionospheric Electrodynamics" ! call IE_Initialize(iError) write(*,*) "Setting nmlts and nlats" ! Aaron Ridley says the original AMIE comes as 1min NH arrays of 25 MLTs (0-24 - NOT 1-25!!!) ! every 2 deg from 90N to 44N (24!)- However, the code interpolates to whatever grid is ! set below as lats and mlts for input with dimensions set in IO_SetnMLT(LAT)s ! Can find the dimensions of AMIE in readAMIEoutput.f90 ! call IO_SetnMLTs(lonmx+1) ! call IO_SetnLats(ithtrns*2+2) write(*,*) "Setting grid" ! call IO_SetGrid(mlts,lats, iError) endif ! if (iamie>=3) then end subroutine init_amie !----------------------------------------------------------------------- subroutine getamie(iyear,iday,iutsec,iprint,istep) ! ! 04/13 bae: ASTRA AMIE is every 5 min (interp), U MI AMIE is every 1 min (no interp in time) ! ! Read AMIE outputs from amie_ncfile file, returning electric potential, ! auroral mean energy and energy flux at current date and time, ! and the data is linearly interpolated to the model time ! gl - 12/07/2002 ! ! use ModKind ! 6/12: for AMIE U MI use hist_module,only: modeltime use input_module,only: amiesh,amienh use cons_module,only: pi,ylonm,ylatm use magfield_module,only: im,jm,dim,djm use dynamo_module,only: mag2geo,nmlat0,phihm use input_module,only: ! from user input | ctpoten, ! cross-cap potential (kV) (e.g., 45.) | power ! hemispheric power (GW) (e.g., 16.) ! ! Args: integer,intent(in) :: iyear,iday,iutsec,iprint,istep ! ! Local: real :: ! | potm(lonp1,jmxm),efxm(lonp1,jmxm),ekvm(lonp1,jmxm), % Resulted in segfault for potm(i,j)!? ! | alat(jmxm),alon(lonp1),alatm(jmxm),alonm(lonp1) | potm(lonmx+1,jmxm),efxm(lonmx+1,jmxm),ekvm(lonmx+1,jmxm), | alat(jmxm),alon(lonmx+1),alatm(jmxm),alonm(lonmx+1) integer :: ier,lw,liw,isign,intpol(2) integer,allocatable :: iw(:) real,allocatable :: w(:) integer :: i,j,ndmon_nl(13),ndmon_lp(13),nonleap integer :: nn, iset, iset1, m, mp1, n, amie_debug,amieday real :: model_ut, denoma, f1, f2, fmin, fmax, maxefx ! real :: del,xmlt,dmlat,dlatm,dlonm,dmltm,rot,dtr ! 6/12: added by Emery for Ridley binary AMIE files integer, dimension(7) :: itime ! real (Real8_) :: rtime ! for U MI files !! real :: rtime integer :: iError, imlt,ilat,jl real :: hpi ! real :: cps,cpn,latnx,latnn,latsx,latsn,potsx,potsn,potnx,potnn ! real :: mltnx,mltnn,mltsx,mltsn integer :: nt,ihem,add60 real :: frac1,frac2 integer :: ip,ipp real :: polekv,polekev,polemwpm2 ! ! Data: GSWM data given at the middle of each month ! -> day of year for this data assuming non-leap year ! J F M A M J J A S O N D J data ndmon_nl ! non_leap | /31,59,90,120,151,181,212,243,273,304,334,365,396/ data ndmon_lp ! leap year | /31,60,91,121,152,182,213,244,274,305,335,366,397/ ! ! if (iprint > 0) then write(6,"(/,72('-'))") write(6,"('GETAMIE:')") write(6,"('Modeltime(1:4)= ',4i5)") modeltime(1:4) write(6,"('Initial requested iyear=',i4,' iday=',i3,' iutsec=', | i10)") iyear,iday,iutsec endif if (istep<=1) then call init_amie endif ! if (istep<=1) then maxefx=50. ! ! Check for leap year nonleap = int(mod(real(iyear),4.)) ! nonleap year /= 0 ! if (iyear.eq.2000) nonleap = 1 ! ! Check times if (modeltime(3).lt.0.or.modeltime(3).gt.59) then write(6,"('>>> getamie: bad imin=',3i4)") modeltime(3) stop 'getamie' endif if (modeltime(4).lt.0.or.modeltime(4).gt.59) then write(6,"('>>> getamie: bad isec=',3i4)") modeltime(4) stop 'getamie' endif ! model_ut = (modeltime(1)-amieday)*24. + float(iutsec)/3600. model_ut = float(iutsec)/3600. if (iamie==2) then ! 04/13 bae: ASTRA AMIE (modeltime=daynum,hr,min,sec) do nt=2,mxtimes add60 = 0 if (amie_dhm(nt-1,2) .ne. amie_dhm(nt,2)) add60 = 60 if (amie_dhm(nt-1,1)==modeltime(1) | .and. amie_dhm(nt-1,2)==modeltime(2) | .and. amie_dhm(nt-1,3) .le. modeltime(3) .and. | amie_dhm(nt,3)+add60 .ge. modeltime(3)) go to 100 if (amie_dhm(nt,1)==modeltime(1) | .and. amie_dhm(nt,2)==modeltime(2) | .and. amie_dhm(nt-1,3)-add60 .le. modeltime(3) .and. | amie_dhm(nt,3) .ge. modeltime(3)) go to 101 enddo write (6,"(1x,'Cannot find modeltime in amie_dhm; modeltime ', | 'amie_dhm(1,mxtimes) =',4i5,1x,3i4,1x,4i4)") modeltime, | (amie_dhm(1,i),i=1,3),mxtimes,(amie_dhm(mxtimes,j),j=1,3) stop ! Use nt and nt-1 for interpolation over 5 min 100 frac1 = (amie_dhm(nt,3)+add60-modeltime(3))/5. go to 102 101 frac1 = (amie_dhm(nt,3)-modeltime(3))/5. 102 frac2 = 1. - frac1 if (frac1<0. .or. frac1>1.) then write (6,"(1x,'frac1,2 bad; modeltime amie_dhm(nt-1,nt) =', | 2f7.2,4i5,1x,4i4,1x,4i4)") frac1,frac2,modeltime,nt-1, | (amie_dhm(nt-1,i),i=1,3),nt,(amie_dhm(nt,j),j=1,3) stop endif if (iprint > 0) then write(6,"('Initial frac1,2,nt-1,nt,amie_dhm =',2f7.2,2i5,6i4)") | frac1,frac2,nt-1,nt,(amie_dhm(nt-1,i),i=1,3), | (amie_dhm(nt,j),j=1,3) endif ! ASTRA AMIE arrays (kV,mw/m2,keV) are read over 0000to2320 mlts from 90to40, -40to-90 mlat ! but are stored in init_amie in order SH -90 to ~-40 mlat, NH +90 to ~+40mlat ! Store in U MI AMIE arrays (V,mw/m2,keV) in order SH -90 to ~-40 mlat, NH +90 to ~+40mlat do i=1,lonmx+1 do j=1,latp1*2 TempPotential(i,j) = frac1*epkv(nt-1,j,i)*1.e+3 | + frac2*epkv(nt,j,i)*1.e+3 AveEOut(i,j) = frac1*ekev(nt-1,j,i) | + frac2*ekev(nt,j,i) EFluxOut(i,j) = frac1*wpm2(nt-1,j,i) | + frac2*wpm2(nt,j,i) enddo enddo endif ! if (iamie==2) then if (iamie>=3) then ! 6/12: added by Emery for Ridley binary AMIE ! modeltime(1) = daynumber of year itime(1) = iyear do i=1,12 if (nonleap == 0) then if (ndmon_lp(i) <= modeltime(1)) then itime(2) = i+1 endif else if (ndmon_nl(i) <= modeltime(1)) then itime(2) = i+1 endif endif enddo if (nonleap == 0) then itime(3) = modeltime(1)-ndmon_lp(itime(2)-1) else itime(3) = modeltime(1)-ndmon_nl(itime(2)-1) endif itime(4) = modeltime(2) itime(5) = modeltime(3) itime(6) = modeltime(4) itime(7) = 0 ! Calling time conversion for Ridley AMIE !! write(*,*) nonleap,itime ! call time_int_to_real(itime, rtime) !! write(*,*) "rtime =",rtime ! call IO_SetTime(rtime) !! write(*,*) "after SetTime rtime =",rtime ! 8/13 bae: Changed in timing_IO_UMI_AMIE.F Method from IE_Closest_ to IE_Interpolate_ for Eflux and AveE ! call get_AMIE_values(rtime) !! write(*,*) "after get_AMIE rtime =",rtime ! ! AveEOut in keV (from ~1-9 keV) ! call IO_GetAveE(AveEOut, iError) ! if (iError /= 0) then ! write(*,*) "Error : ",iError ! else !! write(*,*) "AveEOut itime=",itime !! write(*,*) AveEOut ! endif ! EfluxOut in mW/m2? (~24 mW/m2 max 06349 0030UT - has to be mW/m2!) ! 6/12: Emery - U MI AMIE EFluxOut is in erg/cm2-s or mW/m2 (NOT in W/m2!) ! call IO_GetEFlux(EfluxOut, iError) ! if (iError /= 0) then ! write(*,*) "Error : ",iError ! else !! write(*,*) "EfluxOut itime=",itime !! write(*,*) EfluxOut ! endif ! call IO_GetPotential(TempPotential, iError) ! if (iError /= 0) then ! write(*,*) "Error : ",iError ! stop ! else ! endif ! for if read potential w error ! if (itime(4).eq.0 .and. itime(5).eq.2 .and. itime(3).eq.13) then ! write(*,*) "i,rtime=",itime,rtime ! mlts = 10x 1-25 (same as 25 IEi_HavenMLTS) ! write(*,*) "mlts=",mlts ! lats = 25x 86,82,78,74,70,66,62,58,54,50 (10 mlats - diff from 39 IEi_HavenLats - expected 88,87-50) ! write(*,*) "lats=",lats ! TempPotential in Volts ! For SWMF, have zeros at poles and at i=1 for each j ! write(*,*) "SH pot in V, -90 to -42 mlat" ! write(*,*) ((TempPotential(i,j),i=1,lonmx),j=1,latp1) ! write(*,*) "NH pot in V, 90 to 42 mlat" ! write(*,*) ((TempPotential(i,j),i=1,lonmx),j=latp1+1,latp1*2) ! endif if (iamie==4) then ! 8/13 bae found SWMF 1=181 constant at all mlats from 2-180, so had to be revised as ave of #2+180! do j=1,latp1+latp1 TempPotential(1,j) = 0.5*(TempPotential(2,j) + | TempPotential(lonmx,j)) AveEOut(1,j) = 0.5*(AveEOut(2,j)+AveEOut(lonmx,j)) EfluxOut(1,j) = 0.5*(EfluxOut(2,j)+EfluxOut(lonmx,j)) TempPotential(lonmx+1,j) = TempPotential(1,j) AveEOut(lonmx+1,j) = AveEOut(1,j) EfluxOut(lonmx+1,j) = EfluxOut(1,j) enddo ! TEMP printout - now looks good (except poles still 0!) ! do ipp=1,latp1 ! write (6,"(1x,'i=1,2 pot,kev,eflx at j=',i3,6e12.4)") ipp, ! | TempPotential(1,ipp),AveEOut(1,ipp),EfluxOut(1,ipp), ! | TempPotential(2,ipp),AveEOut(2,ipp),EfluxOut(2,ipp) ! enddo ! Pole values are trash for SWMF, so compute values at j=1 and j=latp1+1 ! Find pole value as average of closest latitude for SWMF (and for AMIE also) ip = 1 ipp = 2 polekv = 0. polekev = 0. polemwpm2 = 0. do i=1,lonmx polekv = polekv + TempPotential(i,ipp) polekev = polekev + AveEOut(i,ipp) polemwpm2 = polemwpm2 + EfluxOut(i,ipp) enddo polekv = polekv/lonmx polekev = polekev/lonmx polemwpm2 = polemwpm2/lonmx do i=1,lonmx+1 TempPotential(i,ip) = polekv AveEOut(i,ip) = polekev EfluxOut(i,ip) = polemwpm2 enddo ip = latp1+1 ipp = latp1+2 polekv = 0. polekev = 0. polemwpm2 = 0. do i=1,lonmx polekv = polekv + TempPotential(i,ipp) polekev = polekev + AveEOut(i,ipp) polemwpm2 = polemwpm2 + EfluxOut(i,ipp) enddo polekv = polekv/lonmx polekev = polekev/lonmx polemwpm2 = polemwpm2/lonmx do i=1,lonmx+1 TempPotential(i,ip) = polekv AveEOut(i,ip) = polekev EfluxOut(i,ip) = polemwpm2 enddo if (istep<3) write (6,"(1x,'polekv,kev,mwpm2 =',3e12.4)") | polekv,polekev,polemwpm2 endif ! if (iamie==4) then endif ! if (iamie>=3) then ! get rid of outragous flux values > maxefx ! Calculate hpi_s,nh_amie in GW do j=1,2 jl = (j-1)*latp1 hpi = 0. do imlt=1,lonmx do ilat=1,latp1 ! if (iprint>0) write(6,"(1x,'jl,imlt,ilat=',3i3)") jl,imlt,ilat if (jl < 1) then hpi = hpi + da(ilat)*EfluxOut(imlt,ilat+jl) else hpi = hpi + da(ilat)*EfluxOut(imlt,ilat+jl) endif enddo enddo if (jl < 1) then hpi_sh_amie = hpi*1.e-12 else ! Convert from mW/m2 * m2(da) * 1.e-12GW/mW hpi_nh_amie = hpi*1.e-12 endif enddo ! for j=1,2 and jl ! Recalculate power for aurora.F power = 0.5*(hpi_sh_amie+hpi_nh_amie) write(6,"(1x,'istep UTh HPS,N,av= ',i9,f8.2,3f7.2)")istep, | model_ut,hpi_sh_amie,hpi_nh_amie,power ! ! gl - 14/07/2002 changed from Emery's original amieterp.f ! ! Change mag coords so latitude extends to the equator in into the opposite ! pole, and so MLT is converted to apex longitues between -180 and +180. ! Put EPOTDUP, EMKEV and WPM2 into S. Hem and N. Hem alike in EPOTM, EKEVM, ! and EFLXM, dimensioned (IMXMP,JMNH) where MLAT goes from -90 to ~-2 ! for ihsoln=1, and from +90 to ~+2 for ihsoln=2 since GRDSET is for ! -90 to the equator. (Will work for +90 to the equator also.) ! AMIE grid goes from -90 to ~-40 for SH and +90 to ~+40 for NH (6/12: w dmlat=1.67deg)) ! TGCM grid goes from -87.5 to -2.5 for SH and +2.5 to +87.5 for NH (6/12: w dglat=5deg) ! Hence, at very end, need to reverse order of NH in whole globe arrays. ! mlongitude starts from 180 degree rot = 250./15. - iutsec/3600. dmltm = 24./float(lonmx) do i=1,lonp1 xmlt = (i-1)*dmltm - rot + 24. xmlt = AMOD(xmlt,24.) m = IFIX(xmlt/dmltm + 1.01) mp1 = m + 1 IF (mp1 > lonp1) mp1 = 2 del = xmlt - (m-1)*dmltm C Initialize arrays around equator do j=latp1+1,ithmx potm(i,j) = 0. potm(i,jmxm+1-j) = 0. ekvm(i,j) = (1.-del)*AveEOut(m,latp1) + | del*AveEOut(mp1,latp1) ekvm(i,jmxm+1-j) = (1.-del)*AveEOut(m,latp1+latp1) + | del*AveEOut(mp1,latp1+latp1) efxm(i,j) = 0. efxm(i,jmxm+1-j) = 0. enddo C Put in AMIE arrays from pole to latp1 do j=1,latp1 potm(i,j) = (1.-del)*TempPotential(m,j) + | del*TempPotential(mp1,j) potm(i,jmxm+1-j) = (1.-del)*TempPotential(m,j+latp1) + | del*TempPotential(mp1,j+latp1) ekvm(i,j) = (1.-del)*AveEOut(m,j) + del*AveEOut(mp1,j) ekvm(i,jmxm+1-j) = (1.-del)*AveEOut(m,j+latp1) + | del*AveEOut(mp1,j+latp1) efxm(i,j) = (1.-del)*EfluxOut(m,j) + del*EfluxOut(mp1,j) efxm(i,jmxm+1-j) = (1.-del)*EfluxOut(m,j+latp1) + | del*EfluxOut(mp1,j+latp1) enddo enddo C Set up coeffs to go between POTM(IMXMP,JMNH) and TIEPOT(IMAXM,JMAXMH)=phihm in dynamo C C **** SET GRID SPACING DLATM, DLONG, DLONM C DMLAT=lat spacing in degrees of AMIE apex grid C 6/12: from amie.nc.out in ~ganglu/amie, have 31 lats w dmlat=1.67 from 90 to 40 mlat C and 37 mlts w dmltm=0.67h from 0 to 24MLT, where lonmx=36 (NOT 37!) C 18*3=54 so 90/54=1.67=dmlat (ithtrns=30 to 40 mlat, ithmx=55 to equator, jmxm=2*ithmx-1) ! dtr = pi/180. ! DMLAT = 180. / FLOAT(jmxm-1) ! DLATM = DMLAT*dtr ! DLONM = 2.*PI/FLOAT(lonmx) ! DMLTM = 24./FLOAT(lonmx) C **** C **** SET ARRAY YLATM (LATITUDE VALUES FOR GEOMAGNETIC GRID C **** alatm(1) = -PI/2. alat(1) = -90. alatm(jmxm) = PI/2. alat(jmxm) = 90. do i = 2,ithmx alat(i) = alat(i-1)+dlatm/dtr alat(jmxm+1-i) = alat(jmxm+2-i)-dlatm/dtr alatm(i) = alatm(i-1)+dlatm alatm(jmxm+1-i) = alatm(jmxm+2-i)-dlatm enddo alon(1) = -pi/dtr alonm(1) = -pi do i=2,lonp1 alon(i) = alon(i-1) + dlonm/dtr alonm(i) = alonm(i-1) + dlonm enddo ! ylatm and ylonm are arrays of latitudes and longitudes of the ! distored magnetic grids in radian - from consdyn.h ! Convert from apex magnetic grid to distorted magnetic grid ! ! Allocate workspace for regrid routine rgrd2.f: lw = nmlonp1+nmlat+2*nmlonp1 if (.not. allocated(w)) allocate(w(lw),stat=ier) if (ier /= 0) write(6,"('>>> horizontal_interp: error allocating', | ' w(lw): lw=',i6,' ier=',i4)") lw,ier liw = nmlonp1 + nmlat if (.not. allocated(iw)) allocate(iw(liw),stat=ier) if (ier /= 0) write(6,"('>>> horzontal_interp: error allocating', | ' iw(liw): liw=',i6,' ier=',i4)") liw,ier intpol(:) = 1 ! linear (not cubic) interp in both dimensions if (alatm(1) > ylatm(1)) alatm(1) = ylatm(1) if (alatm(jmxm) < ylatm(nmlat)) alatm(jmxm) = ylatm(nmlat) if (alonm(1) > ylonm(1)) alonm(1) = ylonm(1) if (alonm(lonp1) < ylonm(nmlonp1)) alonm(lonp1) = ylonm(nmlonp1) call rgrd2(lonp1,jmxm,alonm,alatm,potm,nmlonp1,nmlat, | ylonm,ylatm,phihm,intpol,w,lw,iw,liw,ier) call rgrd2(lonp1,jmxm,alonm,alatm,ekvm,nmlonp1,nmlat, | ylonm,ylatm,tieekv,intpol,w,lw,iw,liw,ier) call rgrd2(lonp1,jmxm,alonm,alatm,efxm,nmlonp1,nmlat, | ylonm,ylatm,tieefx,intpol,w,lw,iw,liw,ier) C Convert from TGCM distorted magnetic grid to geographic one call mag2geo(tieekv(1,1),gekev(1,0),im(1,0),jm(1,0), | dim(1,0),djm(1,0),nlonp1,nmlonp1,nlon,nlat+2,nmlon,nmlat) call mag2geo(tieefx(1,1),gmwpm2(1,0),im(1,0),jm(1,0), | dim(1,0),djm(1,0),nlonp1,nmlonp1,nlon,nlat+2,nmlon,nmlat) C **** C **** INSERT PERIODIC POINTS C **** C 4/13 bae: use i=1,nlon of gekev and gmwpm2 in aurora.F and below, BUT nlonp1 is needed for alfa,eflux! DO j = 1,nlat gekev(nlonp1,j) = gekev(1,j) gmwpm2(nlonp1,j) = gmwpm2(1,j) ENDDO ! temp - save the AMIE output in magnetic grid amie_debug = 1 if (amie_debug > 0) then ! call addfld('TIEPOT',' ',' ',phihm,'mlon',1,nmlonp1, ! | 'mlat',1,nmlat,0) call addfld('TIEEKV',' ',' ',tieekv,'mlon',1,nmlonp1, | 'mlat',1,nmlat,0) call addfld('TIEEFX',' ',' ',tieefx,'mlon',1,nmlonp1, | 'mlat',1,nmlat,0) endif ! if (iprint > 0) then write(6,"('getamie: AMIE data interpolated to date and time')") write(6,"(72('-'),/)") endif end subroutine getamie !------------------------------------------------------------------- end module amie_module