ameriflux/12.plot.ncl
author Forrest Hoffman <forrest@climatemodeling.org>
Mon, 26 Jan 2009 22:08:20 -0500
changeset 0 0c6405ab2ff4
permissions -rw-r--r--
Initial commit of C-LAMP Diagnostics from Jeff Lee
     1 ;************************************************************
     2 ; required command line input parameters:
     3 ;  ncl 'model_name="10cn" model_grid="T42" dirm="/.../ film="..."' 01.npp.ncl
     4 ;
     5 ; using gsn_table for all
     6 ; output: line plot for each site (4 fields)
     7 ;         table for M_score
     8 ;************************************************************
     9 load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/gsn_code.ncl"
    10 load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/gsn_csm.ncl"
    11 load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/contributed.ncl"
    12 ;************************************************************
    13 procedure set_line(lines:string,nline:integer,newlines:string) 
    14 begin
    15 ; add line to ascci/html file
    16     
    17   nnewlines = dimsizes(newlines)
    18   if(nline+nnewlines-1.ge.dimsizes(lines))
    19     print("set_line: bad index, not setting anything.") 
    20     return
    21   end if 
    22   lines(nline:nline+nnewlines-1) = newlines
    23 ;  print ("lines = " + lines(nline:nline+nnewlines-1))
    24   nline = nline + nnewlines
    25   return 
    26 end
    27 ;*************************************************************
    28 begin
    29 
    30   plot_type     = "ps"
    31   plot_type_new = "png"
    32 
    33 ; for 6 fields, 12-monthly
    34   nmon      = 12
    35   nfield    = 6
    36 
    37 ;************************************************
    38 ; read model data
    39 ;************************************************
    40   model = "cn"
    41 ; model = "casa"
    42   
    43   model_name = "i01.10" + model
    44 
    45   dirm  = "/fis/cgd/cseg/people/jeff/surface_data/"
    46   film  = "lnd_T42.nc"
    47   fm    = addfile(dirm+film,"r")
    48   
    49   xm    = fm->lon
    50   ym    = fm->lat
    51 ;------------------------------------------------
    52   nlat = dimsizes(ym)
    53   nlon = dimsizes(xm)
    54 
    55   data_mod0 = new ((/nfield,nmon,nlat,nlon/),float)
    56 
    57 ; change to unit of observed (u mol/m2/s)
    58 ; Model_units [=] gC/m2/s
    59 ; 12. = molecular weight of C
    60 ; u mol = 1e-6 mol
    61   factor = 1e6 /12.
    62 
    63   ENERGY ="new"
    64 
    65 ;************************************************
    66 ; read data: observed
    67 ;************************************************
    68 
    69  station = (/"ARM_Oklahoma" \
    70             ,"ARM_Oklahoma_burn" \
    71             ,"ARM_Oklahoma_control" \
    72             ,"Atqasuk" \
    73             ,"Audubon" \
    74             ,"AustinCary" \
    75             ,"Bartlett" \
    76             ,"Bondville" \
    77             ,"Brookings" \
    78             ,"Donaldson" \
    79             ,"Duke_Forest_Hardwoods" \
    80             ,"Duke_Forest_Open_Field" \
    81             ,"Duke_Forest_Pine" \
    82             ,"Fermi_Ag" \
    83             ,"Fermi_Prairie" \
    84             ,"Flagstaff_Managed" \
    85             ,"Flagstaff_Unmanaged" \
    86             ,"Flagstaff_Wildfire" \
    87             ,"FortPeck" \
    88             ,"FreemanRanch_mesquite" \
    89             ,"Goodwin_Creek" \
    90             ,"HarvardForest" \
    91             ,"HarvardForestHemlock" \
    92             ,"HowlandForestMain" \
    93             ,"HowlandForestWest" \
    94             ,"Ivotuk" \
    95             ,"KendallGrasslands" \
    96             ,"KennedySpaceCenterPine" \
    97             ,"KennedySpaceCenterScrub" \
    98             ,"LittleProspect" \
    99             ,"LostCreek" \
   100             ,"Mead-irrigated" \
   101             ,"Mead-irrigated-rotation" \
   102             ,"Mead-rainfed" \
   103             ,"Metolius_2nd_YoungPonderosaPine" \
   104             ,"MetoliusEyerly" \
   105             ,"MetoliusIntermediatePine" \
   106             ,"MetoliusOldPonderosaPine" \
   107             ,"MissouriOzark" \
   108             ,"Mize" \
   109             ,"MorganMonroe" \
   110             ,"NiwotRidge" \
   111             ,"NorthCarolina_cc" \
   112             ,"NorthCarolina_lp" \
   113             ,"ParkFalls" \
   114             ,"Rayonier" \
   115             ,"SantaRita" \
   116             ,"SkyOaks_Old" \
   117             ,"SkyOaks_PostFire" \
   118             ,"SkyOaks_Young" \
   119             ,"SylvaniaWilderness" \
   120             ,"Toledo" \
   121             ,"Tonzi" \
   122             ,"UCI_1850" \
   123             ,"UCI_1930" \
   124             ,"UCI_1964" \
   125             ,"UCI_1964wet" \
   126             ,"UCI_1981" \
   127             ,"UCI_1989" \
   128             ,"UCI_1998" \
   129             ,"UMBS" \
   130             ,"Vaira" \
   131             ,"WalkerBranch" \
   132             ,"WillowCreek" \
   133             ,"WindRiver" \
   134             ,"Wisconsin_ihw" \
   135             ,"Wisconsin_irp" \
   136             ,"Wisconsin_mrp" \
   137             ,"Wisconsin_myjp" \
   138             ,"Wisconsin_pb" \
   139             ,"Wisconsin_rpcc" \
   140             ,"Wisconsin_yhw" \
   141             ,"Wisconsin_yjp" \
   142             ,"Wisconsin_yrp" \
   143             /)
   144 
   145  year_ob = (/"2003-2006" \
   146             ,"2005-2006" \
   147             ,"2005-2006" \
   148             ,"1999-2006" \
   149             ,"2002-2006" \
   150             ,"2000-2005" \
   151             ,"2004-2005" \
   152             ,"1996-2006" \
   153             ,"2004-2006" \
   154             ,"1999-2004" \
   155             ,"2003-2005" \
   156             ,"2001-2005" \
   157             ,"2001-2005" \
   158             ,"2005-2006" \
   159             ,"2004-2006" \
   160             ,"2005-2006" \
   161             ,"2005-2006" \
   162             ,"2005-2006" \
   163             ,"2000-2006" \
   164             ,"2004-2006" \
   165             ,"2002-2006" \
   166             ,"1991-2004" \
   167             ,"2004-2004" \
   168             ,"1996-2004" \
   169             ,"1999-2004" \
   170             ,"2003-2006" \
   171             ,"2004-2006" \
   172             ,"2002-2002" \
   173             ,"2000-2006" \
   174             ,"2002-2005" \
   175             ,"2001-2005" \
   176             ,"2001-2005" \
   177             ,"2001-2005" \
   178             ,"2001-2005" \
   179             ,"2004-2005" \
   180             ,"2004-2005" \
   181             ,"2003-2005" \
   182             ,"1996-2000" \
   183             ,"2004-2006" \
   184             ,"1998-2004" \
   185             ,"1999-2005" \
   186             ,"1999-2003" \
   187             ,"2005-2006" \
   188             ,"2005-2006" \
   189             ,"1996-2003" \
   190             ,"1998-1998" \
   191             ,"2004-2006" \
   192             ,"1997-2006" \
   193             ,"2004-2006" \
   194             ,"1997-2006" \
   195             ,"2002-2006" \
   196             ,"2004-2005" \
   197             ,"2001-2006" \
   198             ,"2004-2005" \
   199             ,"2001-2005" \
   200             ,"2001-2005" \
   201             ,"2002-2004" \
   202             ,"2001-2005" \
   203             ,"2001-2005" \
   204             ,"2002-2005" \
   205             ,"1999-2003" \
   206             ,"2001-2006" \
   207             ,"1995-1999" \
   208             ,"1999-2005" \
   209             ,"1999-2004" \
   210             ,"2003-2003" \
   211             ,"2003-2003" \
   212             ,"2002-2005" \
   213             ,"2004-2004" \
   214             ,"2002-2002" \
   215             ,"2005-2005" \
   216             ,"2002-2002" \
   217             ,"2004-2005" \
   218             ,"2002-2002" \
   219             /)
   220 
   221  field   = (/"CO2 Flux" \
   222             ,"Net Radiation" \
   223             ,"Latent Heat" \
   224             ,"Sensible Heat" \
   225             ,"GPP Flux" \
   226             ,"Respiration" \
   227             /)
   228 
   229  nstation  = dimsizes(station)
   230 
   231  data_mod  = new ((/nstation, nfield, nmon/),float)
   232  data_ob   = new ((/nstation, nfield, nmon/),float)
   233  lat_ob    = new ((/nstation/),float)
   234  lon_ob    = new ((/nstation/),float)
   235 
   236  diro_root  = "/fis/cgd/cseg/people/jeff/clamp/ameriflux/data/"
   237  dirm_root  = "/fis/cgd/cseg/people/jeff/clamp_data/model/"
   238 
   239  do n = 0,nstation-1
   240 
   241 ;-------------------------------------------------
   242 ;   get ob data
   243 
   244     diro = diro_root + station(n)+"/"
   245     filo = year_ob(n)+"_L4_m.nc"
   246     fo   = addfile (diro+filo,"r")
   247 
   248     print (filo)
   249  
   250     lon_ob(n) = fo->lon 
   251     lat_ob(n) = fo->lat
   252 
   253     data      = fo->NEE_or_fMDS
   254     data_ob(n,0,:) = dim_avg(data(month|:,year|:))
   255     delete (data)
   256 
   257     data      = fo->Rg_f
   258     data_ob(n,1,:) = dim_avg(data(month|:,year|:))
   259     delete (data)
   260 
   261     data      = fo->LE_f
   262     data_ob(n,2,:) = dim_avg(data(month|:,year|:))
   263     delete (data)
   264 
   265     data      = fo->H_f
   266     data_ob(n,3,:) = dim_avg(data(month|:,year|:))
   267     delete (data)
   268 
   269     data      = fo->GPP_or_MDS
   270     data_ob(n,4,:) = dim_avg(data(month|:,year|:))
   271     delete (data)
   272 
   273     data      = fo->Reco_or
   274     data_ob(n,5,:) = dim_avg(data(month|:,year|:))
   275     delete (data)
   276 
   277     delete (fo)
   278 ;---------------------------------------------------
   279 ;   get model data
   280 
   281 ;   film = model_name+"_"+ year_ob(n)+"_MONS_climo.nc"
   282     film = model_name+"_"+"1990-2004"+"_MONS_climo.nc"
   283     fm   = addfile (dirm_root+film,"r")
   284 
   285     print (film)
   286 
   287 if (ENERGY .eq. "old") then
   288 
   289   data = fm->NEE
   290   data_mod0(0,:,:,:) = data(:,:,:) * factor
   291   delete (data)
   292 
   293 ; data  = fm->LATENT
   294   data1 = fm->FCEV
   295   data2 = fm->FCTR
   296   data3 = fm->FGEV
   297   data_mod0(2,:,:,:) = data1(:,:,:)+data2(:,:,:)+data3(:,:,:) 
   298   delete (data1)
   299   delete (data2)
   300   delete (data3)
   301 
   302 ; data = fm->SENSIBLE
   303   data  = fm->FSH
   304   data_mod0(3,:,:,:) = data(:,:,:) 
   305   delete (data)
   306 
   307 ; data  = fm->NETRAD
   308   data1 = fm->FSA
   309   data2 = fm->FIRA
   310   data_mod0(1,:,:,:) = data1(:,:,:)-data2(:,:,:)-data_mod0(2,:,:,:)-data_mod0(3,:,:,:) 
   311   delete (data1)
   312   delete (data2)
   313 
   314 else
   315 
   316   data = fm->NEE
   317   data_mod0(0,:,:,:) = data(:,:,:) * factor
   318   delete (data)
   319 
   320   data = fm->NETRAD
   321   data_mod0(1,:,:,:) = data(:,:,:) 
   322   delete (data)
   323 
   324   data = fm->LATENT
   325   data_mod0(2,:,:,:) = data(:,:,:) 
   326   delete (data)
   327 
   328 ; data = fm->SENSIBLE
   329   data = fm->FSH
   330   data_mod0(3,:,:,:) = data(:,:,:) 
   331   delete (data)
   332 end if
   333 
   334   data = fm->GPP
   335   data_mod0(4,:,:,:) = data(:,:,:) * factor 
   336   delete (data)
   337 
   338   if (model .eq. "cn") then
   339      data = fm->ER
   340   else
   341      data1 = fm->AR
   342      data2 = fm->HR
   343      data  = data1 + data2
   344     
   345      delete (data1)
   346      delete (data2)
   347   end if
   348 
   349   data_mod0(5,:,:,:) = data(:,:,:) * factor
   350   delete (data)
   351 
   352   delete (fm) 
   353 
   354 ;************************************************************
   355 ; interpolate model data into observed station
   356 ; note: model is 0-360E, 90S-90N
   357 ;************************************************************
   358 
   359 ; to be able to handle observation at (-89.98,-24.80)
   360   ym(0) = -90.  
   361 
   362   yy = linint2_points_Wrap(xm,ym,data_mod0,True,lon_ob(n),lat_ob(n),0)
   363 
   364 ; print (yy)
   365  
   366   data_mod(n,:,:) = yy(:,:,0)
   367 
   368   delete (yy)
   369 
   370  end do
   371 ;************************************************************
   372 ; compute correlation coef and M score
   373 ;************************************************************
   374 
   375  score_max = 5.
   376 
   377  ccr     = new ((/nstation, nfield/),float)
   378  M_score = new ((/nstation, nfield/),float) 
   379 
   380  do n=0,nstation-1
   381  do m=0,nfield-1   
   382     ccr(n,m) = esccr(data_ob(n,m,:),data_mod(n,m,:),0)
   383     bias = sum(abs(data_mod(n,m,:)-data_ob(n,m,:))/(abs(data_mod(n,m,:))+abs(data_ob(n,m,:))))
   384     M_score(n,m) = (1. -(bias/nmon)) * score_max
   385  end do
   386  end do
   387 
   388  M_co2 = avg(M_score(:,0))
   389  M_rad = avg(M_score(:,1))
   390  M_lh  = avg(M_score(:,2))
   391  M_sh  = avg(M_score(:,3))
   392  M_gpp = avg(M_score(:,4))
   393  M_er  = avg(M_score(:,5))
   394  M_all = M_co2+ M_rad +M_lh + M_sh + M_gpp + M_er
   395 
   396  M_energy_co2 = sprintf("%.2f", M_co2)
   397  M_energy_rad = sprintf("%.2f", M_rad)
   398  M_energy_lh  = sprintf("%.2f", M_lh )
   399  M_energy_sh  = sprintf("%.2f", M_sh )
   400  M_energy_gpp = sprintf("%.2f", M_gpp)
   401  M_energy_er  = sprintf("%.2f", M_er )
   402  M_energy_all = sprintf("%.2f", M_all)
   403 
   404 ;*******************************************************************
   405 ; for station line plot
   406 ;*******************************************************************
   407 
   408 ; for x-axis in xyplot
   409   mon = ispan(1,12,1)
   410   mon@long_name = "month"
   411 
   412   res                   = True               ; plot mods desired
   413   res@xyLineThicknesses = (/2.0,2.0/)        ; make 2nd lines thicker
   414   res@xyLineColors      = (/"blue","red"/)   ; line color (ob,model)
   415 ;-------------------------------------------------------------------------
   416 ; Add a boxed legend using the more simple method
   417 
   418   res@pmLegendDisplayMode    = "Always"
   419 ; res@pmLegendWidthF         = 0.1
   420   res@pmLegendWidthF         = 0.08
   421   res@pmLegendHeightF        = 0.06
   422 ; res@pmLegendOrthogonalPosF = -1.17
   423 ; res@pmLegendOrthogonalPosF = -1.00  ;(downward)
   424   res@pmLegendOrthogonalPosF = -0.30  ;(downward)
   425 
   426 ; res@pmLegendParallelPosF   =  0.18
   427   res@pmLegendParallelPosF   =  0.23  ;(rightward)
   428 
   429 ; res@lgPerimOn             = False
   430   res@lgLabelFontHeightF     = 0.015
   431   res@xyExplicitLegendLabels = (/"observed",model_name/)
   432 ;-------------------------------------------------------------------
   433 ; for panel plot
   434   res@gsnFrame     = False                   ; Do not draw plot 
   435   res@gsnDraw      = False                   ; Do not advance frame
   436 
   437   pres                            = True     ; panel plot mods desired
   438   pres@gsnPanelYWhiteSpacePercent = 5        ; increase white space around
   439                                              ; indiv. plots in panel
   440   pres@gsnMaximize                = True     ; fill the page
   441 ;-------------------------------------------------------------------
   442 
   443   plot_data   = new((/2,12/),float)
   444   plot_data!0 = "case"
   445   plot_data!1 = "month"
   446 
   447   do n = 0,nstation-1
   448 ;----------------------------
   449 ; for observed
   450 
   451     plot_name = station(n)+"_ob"    
   452     title = station(n)+"("+sprintf("%5.2f",lat_ob(n))+","+sprintf("%5.2f",lon_ob(n))+")"
   453     res@tiMainString = title
   454 
   455     wks = gsn_open_wks (plot_type,plot_name)
   456     plot=new(nfield,graphic)                        ; create graphic array   
   457                            
   458     plot_data(0,:) = (/data_ob (n,0,:)/)
   459     plot_data@long_name = field(0)   
   460     plot(0)=gsn_csm_xy(wks,mon,plot_data(0,:),res)   ; create plot 1
   461 
   462     plot_data(0,:) = (/data_ob (n,1,:)/)
   463     plot_data@long_name = field(1)
   464     plot(1)=gsn_csm_xy(wks,mon,plot_data(0,:),res)   ; create plot 2
   465 
   466     plot_data(0,:) = (/data_ob (n,2,:)/)
   467     plot_data@long_name = field(2)   
   468     plot(2)=gsn_csm_xy(wks,mon,plot_data(0,:),res)   ; create plot 3
   469 
   470     plot_data(0,:) = (/data_ob (n,3,:)/)
   471     plot_data@long_name = field(3)
   472     plot(3)=gsn_csm_xy(wks,mon,plot_data(0,:),res)   ; create plot 4
   473 
   474     plot_data(0,:) = (/data_ob (n,4,:)/)
   475     plot_data@long_name = field(4)
   476     plot(4)=gsn_csm_xy(wks,mon,plot_data(0,:),res)   ; create plot 5
   477 
   478     plot_data(0,:) = (/data_ob (n,5,:)/)
   479     plot_data@long_name = field(5)
   480     plot(5)=gsn_csm_xy(wks,mon,plot_data(0,:),res)   ; create plot 6
   481 
   482     gsn_panel(wks,plot,(/3,2/),pres)                 ; create panel plot
   483 
   484     system("convert "+plot_name+"."+plot_type+" "+plot_name+"."+plot_type_new+";"+ \
   485            "rm "+plot_name+"."+plot_type)
   486 
   487     clear (wks)  
   488     delete (plot)
   489 ;----------------------------
   490 ; for model_vs_ob
   491 
   492     plot_name = station(n)+"_model_vs_ob"
   493     title = station(n)+"("+sprintf("%5.2f",lat_ob(n))+","+sprintf("%5.2f",lon_ob(n))+")"    
   494     res@tiMainString = title
   495 
   496     wks = gsn_open_wks (plot_type,plot_name)
   497     plot=new(nfield,graphic)                        ; create graphic array   
   498                            
   499     plot_data(0,:) = (/data_ob (n,0,:)/)
   500     plot_data(1,:) = (/data_mod(n,0,:)/)
   501     plot_data@long_name = field(0)   
   502     plot(0)=gsn_csm_xy(wks,mon,plot_data,res)   ; create plot 1
   503 
   504     plot_data(0,:) = (/data_ob (n,1,:)/)
   505     plot_data(1,:) = (/data_mod(n,1,:)/)
   506     plot_data@long_name = field(1)
   507     plot(1)=gsn_csm_xy(wks,mon,plot_data,res)   ; create plot 2
   508 
   509     plot_data(0,:) = (/data_ob (n,2,:)/)
   510     plot_data(1,:) = (/data_mod(n,2,:)/)
   511     plot_data@long_name = field(2)   
   512     plot(2)=gsn_csm_xy(wks,mon,plot_data,res)   ; create plot 3
   513 
   514     plot_data(0,:) = (/data_ob (n,3,:)/)
   515     plot_data(1,:) = (/data_mod(n,3,:)/)
   516     plot_data@long_name = field(3)
   517     plot(3)=gsn_csm_xy(wks,mon,plot_data,res)   ; create plot 4
   518 
   519 
   520     plot_data(0,:) = (/data_ob (n,4,:)/)
   521     plot_data(1,:) = (/data_mod(n,4,:)/)
   522     plot_data@long_name = field(4)
   523     plot(4)=gsn_csm_xy(wks,mon,plot_data,res)   ; create plot 5
   524 
   525 
   526     plot_data(0,:) = (/data_ob (n,5,:)/)
   527     plot_data(1,:) = (/data_mod(n,5,:)/)
   528     plot_data@long_name = field(5)
   529     plot(5)=gsn_csm_xy(wks,mon,plot_data,res)   ; create plot 6
   530 
   531     gsn_panel(wks,plot,(/3,2/),pres)                 ; create panel plot
   532 
   533     system("convert "+plot_name+"."+plot_type+" "+plot_name+"."+plot_type_new+";"+ \
   534            "rm "+plot_name+"."+plot_type)
   535 
   536     clear (wks)  
   537     delete (plot)
   538  end do
   539 
   540 ;*******************************************************************
   541 ; html table of site: observed
   542 ;*******************************************************************
   543   output_html = "line_ob.html"
   544 
   545   header = (/"<HTML>" \
   546             ,"<HEAD>" \
   547             ,"<TITLE>CLAMP metrics</TITLE>" \
   548             ,"</HEAD>" \
   549             ,"<H1>Energy at Site: Observation</H1>" \
   550             /) 
   551   footer = "</HTML>"
   552 
   553   table_header = (/ \
   554         "<table border=1 cellspacing=0 cellpadding=3 width=60%>" \
   555        ,"<tr>" \
   556        ,"   <th bgcolor=DDDDDD >Site Name</th>" \
   557        ,"   <th bgcolor=DDDDDD >Latitude</th>" \
   558        ,"   <th bgcolor=DDDDDD >Longitude</th>" \
   559        ,"   <th bgcolor=DDDDDD >Observed</th>" \ 
   560        ,"</tr>" \
   561        /)
   562   table_footer = "</table>"
   563   row_header = "<tr>"
   564   row_footer = "</tr>"
   565 
   566   lines = new(50000,string)
   567   nline = 0
   568 
   569   set_line(lines,nline,header)
   570   set_line(lines,nline,table_header)
   571 ;-----------------------------------------------
   572 ; row of table
   573   
   574   do n = 0,nstation-1
   575      set_line(lines,nline,row_header)
   576 
   577      txt0 = station(n)
   578      txt1 = sprintf("%5.2f", lat_ob(n))
   579      txt2 = sprintf("%5.2f", lon_ob(n))
   580      txt3 = year_ob(n)
   581 
   582      set_line(lines,nline,"<th><a href="+txt0+"_ob.png>"+txt0+"</a></th>")
   583      set_line(lines,nline,"<th>"+txt1+"</th>")
   584      set_line(lines,nline,"<th>"+txt2+"</th>")
   585      set_line(lines,nline,"<th>"+txt3+"</th>")
   586 
   587      set_line(lines,nline,row_footer)
   588   end do
   589 ;-----------------------------------------------
   590   set_line(lines,nline,table_footer)
   591   set_line(lines,nline,footer) 
   592 
   593 ; Now write to an HTML file.
   594   idx = ind(.not.ismissing(lines))
   595   if(.not.any(ismissing(idx))) then
   596     asciiwrite(output_html,lines(idx))
   597   else
   598    print ("error?")
   599   end if
   600   delete (idx)
   601 
   602 ;*******************************************************************
   603 ; score and line table : model vs observed
   604 ;*******************************************************************
   605   output_html = "score+line_vs_ob.html"
   606 
   607   header = (/"<HTML>" \
   608             ,"<HEAD>" \
   609             ,"<TITLE>CLAMP metrics</TITLE>" \
   610             ,"</HEAD>" \
   611             ,"<H1>Energy at Site: Model "+model_name+"</H1>" \
   612             /) 
   613   footer = "</HTML>"
   614 
   615   delete (table_header)
   616   table_header = (/ \
   617         "<table border=1 cellspacing=0 cellpadding=3 width=100%>" \
   618        ,"<tr>" \
   619        ,"   <th bgcolor=DDDDDD >Site Name</th>" \
   620        ,"   <th bgcolor=DDDDDD >Latitude</th>" \
   621        ,"   <th bgcolor=DDDDDD >Longitude</th>" \
   622        ,"   <th bgcolor=DDDDDD >Observed</th>" \
   623        ,"   <th bgcolor=DDDDDD >CO2 Flux</th>" \
   624        ,"   <th bgcolor=DDDDDD >Net Radiation</th>" \
   625        ,"   <th bgcolor=DDDDDD >Latent Heat</th>" \
   626        ,"   <th bgcolor=DDDDDD >Sensible Heat</th>" \
   627        ,"   <th bgcolor=DDDDDD >GPP Glux</th>" \
   628        ,"   <th bgcolor=DDDDDD >Respiration</th>" \
   629        ,"   <th bgcolor=DDDDDD >Average</th>" \
   630        ,"</tr>" \
   631        /)
   632   table_footer = "</table>"
   633   row_header = "<tr>"
   634   row_footer = "</tr>"
   635 
   636   lines = new(50000,string)
   637   nline = 0
   638 
   639   set_line(lines,nline,header)
   640   set_line(lines,nline,table_header)
   641 ;-----------------------------------------------
   642 ; row of table
   643   
   644   do n = 0,nstation-1
   645      set_line(lines,nline,row_header)
   646 
   647      txt0  = station(n)
   648      txt1  = sprintf("%5.2f", lat_ob(n))
   649      txt2  = sprintf("%5.2f", lon_ob(n))
   650      txt3  = year_ob(n)
   651      txt4  = sprintf("%5.2f", M_score(n,0))
   652      txt5  = sprintf("%5.2f", M_score(n,1))
   653      txt6  = sprintf("%5.2f", M_score(n,2))
   654      txt7  = sprintf("%5.2f", M_score(n,3))
   655      txt8  = sprintf("%5.2f", M_score(n,4))
   656      txt9  = sprintf("%5.2f", M_score(n,5))
   657      txt10 = sprintf("%5.2f", avg(M_score(n,:)))
   658 
   659      set_line(lines,nline,"<th><a href="+txt0+"_model_vs_ob.png>"+txt0+"</a></th>")
   660      set_line(lines,nline,"<th>"+txt1+"</th>")
   661      set_line(lines,nline,"<th>"+txt2+"</th>")
   662      set_line(lines,nline,"<th>"+txt3+"</th>")
   663      set_line(lines,nline,"<th>"+txt4+"</th>")
   664      set_line(lines,nline,"<th>"+txt5+"</th>")
   665      set_line(lines,nline,"<th>"+txt6+"</th>")
   666      set_line(lines,nline,"<th>"+txt7+"</th>")
   667      set_line(lines,nline,"<th>"+txt8+"</th>")
   668      set_line(lines,nline,"<th>"+txt9+"</th>")
   669      set_line(lines,nline,"<th>"+txt10+"</th>")
   670 
   671      set_line(lines,nline,row_footer)
   672   end do
   673 
   674 ; last row, summary
   675   set_line(lines,nline,row_header)
   676 
   677   txt0  = "All_"+sprintf("%.0f", nstation)
   678   txt1  = "-"
   679   txt2  = "-"
   680   txt3  = "-"
   681   txt4  = M_energy_co2
   682   txt5  = M_energy_rad
   683   txt6  = M_energy_lh
   684   txt7  = M_energy_sh
   685   txt8  = M_energy_gpp
   686   txt9  = M_energy_er
   687   txt10 = M_energy_all
   688 
   689   set_line(lines,nline,"<th>"+txt0+"</th>")
   690   set_line(lines,nline,"<th>"+txt1+"</th>")
   691   set_line(lines,nline,"<th>"+txt2+"</th>")
   692   set_line(lines,nline,"<th>"+txt3+"</th>")
   693   set_line(lines,nline,"<th>"+txt4+"</th>")
   694   set_line(lines,nline,"<th>"+txt5+"</th>")
   695   set_line(lines,nline,"<th>"+txt6+"</th>")
   696   set_line(lines,nline,"<th>"+txt7+"</th>")
   697   set_line(lines,nline,"<th>"+txt8+"</th>")
   698   set_line(lines,nline,"<th>"+txt9+"</th>")
   699   set_line(lines,nline,"<th>"+txt10+"</th>")
   700 
   701   set_line(lines,nline,row_footer)
   702 ;-----------------------------------------------
   703   set_line(lines,nline,table_footer)
   704   set_line(lines,nline,footer) 
   705 
   706 ; Now write to an HTML file.
   707   idx = ind(.not.ismissing(lines))
   708   if(.not.any(ismissing(idx))) then
   709     asciiwrite(output_html,lines(idx))
   710   else
   711    print ("error?")
   712   end if
   713   delete (idx)
   714 
   715 ;***************************************************************************
   716 ; output plots
   717 ;***************************************************************************
   718   output_dir = model_name+"/ameriflux"
   719 
   720   system("mv *.png *.html " + output_dir) 
   721 ;***************************************************************************
   722 end