PISM, A Parallel Ice Sheet Model  stable v0.5
src/earth/tryLCbd.cc
Go to the documentation of this file.
00001 // Copyright (C) 2007--2011 Ed Bueler and Constantine Khroulev
00002 //
00003 // This file is part of PISM.
00004 //
00005 // PISM is free software; you can redistribute it and/or modify it under the
00006 // terms of the GNU General Public License as published by the Free Software
00007 // Foundation; either version 2 of the License, or (at your option) any later
00008 // version.
00009 //
00010 // PISM is distributed in the hope that it will be useful, but WITHOUT ANY
00011 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00012 // FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00013 // details.
00014 //
00015 // You should have received a copy of the GNU General Public License
00016 // along with PISM; if not, write to the Free Software
00017 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00018 
00019 // get this message with 'tryLCbd -help':
00020 static char help[] =
00021   "Simple testing program for Lingle & Clark bed deformation model.\n"
00022   "  Runs go for 150,000 years on 63.5km grid with 100a time steps and Z=2 in L&C model.\n"
00023   "  SCENARIOS:  run './tryLCbd N' where N=1,2,3,4 as follows\n"
00024   "     (1) dump ice disc on initially level, non-uplifting land, use only viscous \n"
00025   "         half-space model:\n"
00026   "               include_elastic = FALSE, do_uplift = FALSE, H0 = 1000.0\n"
00027   "         center depth b(0,0) should eventually equilibriate to near\n"
00028   "         -1000 * (910/3300) = -275.76 m\n"
00029   "     (2) dump ice disc on initially level, non-uplifting land, use both viscous \n"
00030   "         half-space model and elastic model\n"
00031   "               include_elastic = TRUE, do_uplift = FALSE, H0 = 1000.0\n"
00032   "     (3) never loaded, initially level, uplifting land, use only viscous \n"
00033   "         half-space model (because elastic model gives no additional when no load):\n"
00034   "               include_elastic = FALSE, do_uplift = TRUE, H0 = 0.0\n"
00035   "     (4) dump ice disc on initially level, uplifting land, use both viscous \n"
00036   "         half-space model and elastic model:\n"
00037   "               include_elastic = TRUE, do_uplift = TRUE, H0 = 1000.0\n\n";
00038 
00039 
00040 #include <cmath>
00041 #include <cstdio>
00042 #include <petscvec.h>
00043 #include <petscdmda.h>
00044 #include "pism_const.hh"
00045 #include "NCVariable.hh"
00046 #include "deformation.hh"
00047 #include "pism_options.hh"
00048 
00049 int main(int argc, char *argv[]) {
00050   PetscErrorCode  ierr;
00051 
00052   MPI_Comm    com;  // won't be used except for rank,size
00053   PetscMPIInt rank, size;
00054 
00055   ierr = PetscInitialize(&argc, &argv, PETSC_NULL, help); CHKERRQ(ierr);
00056 
00057   com = PETSC_COMM_WORLD;
00058   ierr = MPI_Comm_rank(com, &rank); CHKERRQ(ierr);
00059   ierr = MPI_Comm_size(com, &size); CHKERRQ(ierr);
00060   
00061   /* This explicit scoping forces destructors to be called before PetscFinalize() */
00062   {
00063     NCConfigVariable config, overrides;
00064 
00065     ierr = init_config(com, rank, config, overrides); CHKERRQ(ierr);
00066 
00067     BedDeformLC bdlc;
00068     DM          da2;
00069     Vec         H, bed, Hstart, bedstart, uplift;
00070     
00071     PetscBool  include_elastic = PETSC_FALSE,
00072                 do_uplift = PETSC_FALSE;
00073     PetscScalar H0 = 1000.0;            // ice disc load thickness
00074 
00075     if (argc >= 2) {
00076       // FIXME:  should use PETSC-style options
00077       switch (argv[1][0]) {
00078         case '1':
00079           include_elastic = PETSC_FALSE;  do_uplift = PETSC_FALSE;  H0 = 1000.0;
00080           break;
00081         case '2':
00082           include_elastic = PETSC_TRUE;  do_uplift = PETSC_FALSE;  H0 = 1000.0;
00083           break;
00084         case '3':
00085           include_elastic = PETSC_FALSE;  do_uplift = PETSC_TRUE;  H0 = 0.0;
00086           break;
00087         case '4':
00088           include_elastic = PETSC_TRUE;  do_uplift = PETSC_TRUE;  H0 = 1000.0;
00089           break;
00090         default:
00091           break; // accept default which is scenario 1
00092       }
00093     }
00094     const PetscScalar R0 = 1000.0e3;          // ice disc load radius
00095     const PetscScalar tfinalyears = 150.0e3;  // total run time
00096 
00097     // FIXME: should accept options here
00098     const PetscInt    Mx = 193, 
00099                       My = 129;
00100 
00101     const PetscScalar Lx = 3000.0e3, 
00102                       Ly = 2000.0e3;
00103     const PetscInt    Z = 2;
00104     const PetscScalar dtyears = 100.0;
00105     
00106     if (rank == 0) { // only runs on proc 0; all sequential
00107       // allocate the variables needed before BedDeformLC can work:
00108       ierr = VecCreateSeq(PETSC_COMM_SELF, Mx*My, &H); CHKERRQ(ierr);
00109       ierr = VecDuplicate(H, &Hstart); CHKERRQ(ierr);
00110       ierr = VecDuplicate(H, &bedstart); CHKERRQ(ierr);
00111       ierr = VecDuplicate(H, &uplift); CHKERRQ(ierr);
00112 
00113       // in order to show bed elevation as a picture, create a da 
00114       ierr = DMDACreate2d(PETSC_COMM_SELF,
00115                           DMDA_BOUNDARY_PERIODIC, DMDA_BOUNDARY_PERIODIC,
00116                           DMDA_STENCIL_STAR,
00117                           My, Mx, PETSC_DECIDE, PETSC_DECIDE, 1, 0,
00118                           PETSC_NULL, PETSC_NULL, &da2); CHKERRQ(ierr);
00119       ierr = DMDASetUniformCoordinates(da2, -Ly, Ly, -Lx, Lx, 
00120                                        PETSC_NULL, PETSC_NULL); CHKERRQ(ierr);
00121       ierr = DMCreateGlobalVector(da2, &bed); CHKERRQ(ierr);
00122 
00123       // create a bed viewer
00124       PetscViewer viewer;
00125       PetscDraw   draw;
00126       const PetscInt  windowx = 500,
00127                       windowy = (PetscInt) (((float) windowx) * Ly / Lx);
00128       ierr = PetscViewerDrawOpen(PETSC_COMM_SELF, PETSC_NULL, "bed elev (m)",
00129            PETSC_DECIDE, PETSC_DECIDE, windowy, windowx, &viewer);  CHKERRQ(ierr);
00130       // following should be redundant, but may put up a title even under 2.3.3-p1:3 where
00131       // there is a no-titles bug
00132       ierr = PetscViewerDrawGetDraw(viewer,0,&draw); CHKERRQ(ierr);
00133       ierr = PetscDrawSetDoubleBuffer(draw); CHKERRQ(ierr);  // remove flicker while we are at it
00134       ierr = PetscDrawSetTitle(draw,"bed elev (m)"); CHKERRQ(ierr);
00135 
00136       // make disc load
00137       ierr = PetscPrintf(PETSC_COMM_SELF,"creating disc load\n"); CHKERRQ(ierr);
00138       // see "Results: Earth deformation only" section of Bueler et al "Fast computation ..."
00139       const PetscScalar dx = (2.0*Lx)/((PetscScalar) Mx - 1), 
00140                         dy = (2.0*Ly)/((PetscScalar) My - 1);
00141       const PetscInt    imid = (Mx-1)/2, jmid = (My-1)/2;
00142       PetscScalar **HH;
00143       ierr = VecGetArray2d(H, Mx, My, 0, 0, &HH); CHKERRQ(ierr);
00144       for (PetscInt i=0; i<Mx; i++) {
00145         for (PetscInt j=0; j<My; j++) {
00146           const PetscScalar r = sqrt( PetscSqr(dx * (i - imid)) + PetscSqr(dy * (j - jmid)) );
00147           if (r < R0) {
00148             HH[i][j] = H0;
00149           } else {
00150             HH[i][j] = 0.0;
00151           }
00152         }
00153       }
00154       ierr = VecRestoreArray2d(H, Mx, My, 0, 0, &HH); CHKERRQ(ierr);
00155       ierr = VecSet(Hstart, 0.0); CHKERRQ(ierr);    // load was zero up till t=0
00156       ierr = VecSet(bedstart, 0.0); CHKERRQ(ierr);       // initially flat bed
00157       
00158       // initialize uplift
00159       if (do_uplift == PETSC_TRUE) {
00160         PetscScalar **upl;
00161         ierr = VecGetArray2d(uplift, Mx, My, 0, 0, &upl); CHKERRQ(ierr);
00162         for (PetscInt i=0; i<Mx; i++) {
00163           for (PetscInt j=0; j<My; j++) {
00164             const PetscScalar peak_up = convert(10, "mm/year", "m/s");  // 10 mm/a
00165             const PetscScalar r = sqrt( PetscSqr(dx * (i - imid)) + PetscSqr(dy * (j - jmid)) );
00166             if (r < 1.5 * R0) {
00167               upl[i][j] = peak_up * (cos(pi * (r / (1.5 * R0))) + 1.0) / 2.0; 
00168             } else {
00169               upl[i][j] = 0.0;
00170             }
00171           }
00172         }
00173         ierr = VecRestoreArray2d(uplift, Mx, My, 0, 0, &upl); CHKERRQ(ierr);
00174       } else {
00175         ierr = VecSet(uplift, 0.0); CHKERRQ(ierr);
00176       }
00177 
00178       ierr = PetscPrintf(PETSC_COMM_SELF,"setting BedDeformLC\n"); CHKERRQ(ierr);
00179       ierr = bdlc.settings(
00180                config,
00181                include_elastic,Mx,My,dx,dy,Z,
00182                &Hstart, &bedstart, &uplift, &H, &bed); CHKERRQ(ierr);
00183 
00184       ierr = PetscPrintf(PETSC_COMM_SELF,"allocating BedDeformLC\n"); CHKERRQ(ierr);
00185       ierr = bdlc.alloc(); CHKERRQ(ierr);
00186       
00187       ierr = PetscPrintf(PETSC_COMM_SELF,"initializing BedDeformLC from uplift map\n"); CHKERRQ(ierr);
00188       ierr = bdlc.uplift_init(); CHKERRQ(ierr);
00189       
00190       ierr = PetscPrintf(PETSC_COMM_SELF,"stepping BedDeformLC\n"); CHKERRQ(ierr);
00191       const PetscInt     KK = (PetscInt) (tfinalyears / dtyears);
00192       PetscScalar **b;
00193       ierr = VecGetArray2d(bedstart, Mx, My, 0, 0, &b); CHKERRQ(ierr);
00194       PetscScalar b0old = b[imid][jmid];
00195       ierr = VecRestoreArray2d(bedstart, Mx, My, 0, 0, &b); CHKERRQ(ierr);
00196       for (PetscInt k=0; k<KK; k++) {
00197         const PetscScalar tyears = k*dtyears;
00198         ierr = bdlc.step(dtyears, tyears); CHKERRQ(ierr);
00199         ierr = VecView(bed,viewer); CHKERRQ(ierr);
00200         ierr = VecGetArray2d(bed, Mx, My, 0, 0, &b); CHKERRQ(ierr);
00201         const PetscScalar b0new = b[imid][jmid];
00202         ierr = VecRestoreArray2d(bed, Mx, My, 0, 0, &b); CHKERRQ(ierr);
00203         const PetscScalar dbdt0 = (b0new - b0old) / (dtyears);
00204 
00205         ierr = PetscPrintf(PETSC_COMM_SELF,
00206                   "   t=%8.0f (a)   b(0,0)=%11.5f (m)  dbdt(0,0)=%11.7f (m/a)\n",
00207                   tyears, b0new, dbdt0); CHKERRQ(ierr);
00208 
00209         char title[100];
00210         snprintf(title,100, "bed elev (m)  [t = %9.1f]", tyears);
00211         ierr = PetscDrawSetTitle(draw,title); CHKERRQ(ierr);
00212         b0old = b0new;
00213       }
00214       ierr = PetscPrintf(PETSC_COMM_SELF,"\ndone\n"); CHKERRQ(ierr);
00215 
00216       ierr = VecDestroy(&H); CHKERRQ(ierr);
00217       ierr = VecDestroy(&bed); CHKERRQ(ierr);
00218       ierr = VecDestroy(&Hstart); CHKERRQ(ierr);
00219       ierr = VecDestroy(&bedstart); CHKERRQ(ierr);
00220       ierr = VecDestroy(&uplift); CHKERRQ(ierr);
00221       ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
00222       ierr = DMDestroy(&da2); CHKERRQ(ierr);
00223     }
00224   }
00225 
00226   ierr = PetscFinalize(); CHKERRQ(ierr);
00227   return 0;
00228 }
 All Classes Namespaces Files Functions Variables Typedefs Friends Defines