Main Page   Class Hierarchy   Compound List   File List   Compound Members   File Members   Related Pages  

jacobi_poisson/JacobiPoissonX.C.sed

Go to the documentation of this file.
00001 /*
00002   Id:           $Id: JacobiPoissonX.C.sed,v 1.1.2.2 2002/03/20 00:43:45 gunney Exp $
00003   Copyright:    (c) 1997-2000 The Regents of the University of California
00004   Release:      $Name: v_0_1_7 $
00005   Revision:     $Revision: 1.1.2.2 $
00006   Modified:     $Date: 2002/03/20 00:43:45 $
00007   Description:  JacobiPoissonX class implementation
00008 */
00009 
00010 
00011 #include <setArrayData.h>
00012 #include <printObjectX.h>
00013 #include <MultiDimArrayAccess.h>
00014 #include <arrayConversionX.h>
00015 #include <setVectorWeightsX.h>
00016 #include "miscFcnsX.h"
00017 #include "JacobiPoissonX.h"
00018 
00019 #include <hier_VariableDatabase2.h>
00020 #include <pdat_CellData2.h>
00021 #include <pdat_SideData2.h>
00022 #include <hier_Patch2.h>
00023 #include <hier_Variable2.h>
00024 #include <geom_CartesianPatchGeometry2.h>
00025 #include <xfer_Geometry2.h>
00026 #include <xfer_CoarsenAlgorithm2.h>
00027 #include <xfer_CoarsenOperator2.h>
00028 #include <xfer_CoarsenSchedule2.h>
00029 #include <xfer_RefineAlgorithm2.h>
00030 #include <xfer_RefineOperator2.h>
00031 #include <xfer_RefineSchedule2.h>
00032 #include <hier_Index2.h>
00033 #include <geom_CartesianCellDoubleLinearRefine2.h>
00034 #include <geom_CartesianCellDoubleConservativeLinearRefine2.h>
00035 #include <geom_CartesianCellDoubleWeightedAverage2.h>
00036 
00037 
00038 using namespace SAMRAI;
00039 
00040 
00041 JacobiPoissonX::JacobiPoissonX(
00042   const string &object_name
00043 , tbox_Pointer<tbox_Database> database
00044 )
00045 : d_context_persistent("PERSISTENT")
00046 , d_context_scratch("SCRATCH")
00047 , d_diffcoef("solution:diffcoef",1)
00048 , d_scalar("solution:scalar",1)
00049 , d_source("source",1)
00050 , d_exact("solution:exact",1)
00051 , d_weight("vector weight",1)
00052 {
00053 
00054   /*
00055     Register variables with VariableDatabase
00056     and get the descriptor indices for those variables.
00057     It is not necessary to save the indices obtained from
00058     the registration, because they can always be retrieved
00059     from the mapVariableAndContextToIndex, but we do it
00060     because we refer to the indices often.
00061   */
00062   {
00063     tbox_Pointer<hier_VariableContext>
00064       context_persistent_ptr(&d_context_persistent,false);
00065     tbox_Pointer<hier_VariableContext>
00066       context_scratch_ptr(&d_context_scratch,false);
00067     hier_VariableDatabaseX *variable_db = hier_VariableDatabaseX::getDatabase();
00068 
00069     /*
00070         Persistent data.
00071     */
00072     d_diffcoef_persistent
00073       = variable_db->registerVariableAndContext(
00074         tbox_Pointer<hier_VariableX>(&d_diffcoef,false)
00075         , context_persistent_ptr );
00076     d_scalar_persistent
00077       = variable_db->registerVariableAndContext(
00078         tbox_Pointer<hier_VariableX>(&d_scalar,false)
00079         , context_persistent_ptr );
00080     d_source_persistent
00081       = variable_db->registerVariableAndContext(
00082         tbox_Pointer<hier_VariableX>(&d_source,false)
00083         , context_persistent_ptr
00084         , hier_IntVectorX(1) /* ghost cell width is 1 */
00085         );
00086     d_exact_persistent
00087       = variable_db->registerVariableAndContext(
00088         tbox_Pointer<hier_VariableX>(&d_exact,false)
00089         , context_persistent_ptr );
00090     d_weight_persistent
00091       = variable_db->registerVariableAndContext(
00092         tbox_Pointer<hier_VariableX>(&d_weight,false)
00093         , context_persistent_ptr );
00094     /*
00095       Scratch data.
00096     */
00097     d_scalar_scratch
00098       = variable_db->registerVariableAndContext(
00099         tbox_Pointer<hier_VariableX>(&d_scalar,false)
00100         , context_scratch_ptr
00101         , hier_IntVectorX(1) /* ghost cell width is 1 */
00102         );
00103   }
00104 
00105   return;
00106 }
00107 
00108 
00109 
00110 
00111 
00112 void JacobiPoissonX::initializeLevelData (
00114     const tbox_Pointer< hier_PatchHierarchyX > hierarchy
00115   , 
00116     const int level_number
00117   , const double init_data_time
00118   , const bool can_be_refined
00119   , 
00120     const bool initial_time
00121   , 
00122     const tbox_Pointer< hier_PatchLevelX > old_level
00123   , const bool allocate_data
00124   ) {
00125 
00126   /*
00127     Reference the level object with the given index from the hierarchy.
00128   */
00129   tbox_Pointer<hier_PatchLevelX> level
00130     = hierarchy->getPatchLevel(level_number);
00131 
00132   /*
00133     If instructed, allocate all patch data on the level.
00134     Allocate only persistent data.  Scratch data will
00135     generally be allocated and deallocated as needed.
00136   */
00137   if ( allocate_data ) {
00138     level->allocatePatchData( d_diffcoef_persistent );
00139     level->allocatePatchData( d_scalar_persistent );
00140     level->allocatePatchData( d_source_persistent );
00141     level->allocatePatchData( d_exact_persistent );
00142     level->allocatePatchData( d_weight_persistent );
00143   }
00144 
00145   /*
00146     Initialize data in all patches in the level.
00147   */
00148   int npatch = level->getNumberOfPatches();
00149   for ( int j=0; j<npatch; j++ ) {
00150 
00151     tbox_Pointer<hier_PatchX> patch
00152       = level->getPatch(j);
00153     tbox_Pointer<geom_CartesianPatchGeometryX> patch_geom
00154       = patch->getPatchGeometry();
00155 
00156     tbox_Pointer<pdat_SideDataX<double> > diffcoef_data
00157       = patch->getPatchData(d_diffcoef_persistent);
00158     tbox_Pointer<pdat_CellDataX<double> > scalar_data
00159       = patch->getPatchData(d_scalar_persistent);
00160     tbox_Pointer<pdat_CellDataX<double> > exact_data
00161       = patch->getPatchData(d_exact_persistent);
00162     tbox_Pointer<pdat_CellDataX<double> > source_data
00163       = patch->getPatchData(d_source_persistent);
00164 
00165     /* Set initial guess of solution. */
00166     MultiDimArrayAccess<double,2> t0
00167       = arrayData2ArrayAccess( scalar_data->getArrayData() );
00168     setArrayDataToConstant(
00169         t0 // arrayData2ArrayAccess( scalar_data->getArrayData() )
00170       , scalar_data->getGhostBox().lower()
00171       , scalar_data->getGhostBox().upper()
00172       , patch_geom->getXLower()
00173       , patch_geom->getXUpper()
00174       , patch_geom->getDx()
00175       , 0.0
00176       );
00177     /* Set exact solution. */
00178     MultiDimArrayAccess<double,2> t4
00179       = arrayData2ArrayAccess( exact_data->getArrayData() );
00180     setArrayDataToSinusoidal(
00181         t4 // arrayData2ArrayAccess( exact_data->getArrayData() )
00182       , exact_data->getBox().lower()
00183       , exact_data->getBox().upper()
00184       , patch_geom->getXLower()
00185       , patch_geom->getXUpper()
00186       , patch_geom->getDx()
00187       );
00188     setArrayDataToScaled(
00189         t4 // arrayData2ArrayAccess( scratch_solution->getArrayData() )
00190       , exact_data->getGhostBox().lower()
00191       , exact_data->getGhostBox().upper()
00192       , -1.0/(8*M_PI*M_PI)
00193       );
00194     /* Set source function. */
00195     MultiDimArrayAccess<double,2> t1
00196       = arrayData2ArrayAccess( source_data->getArrayData() );
00197     setArrayDataToSinusoidal(
00198         t1 // arrayData2ArrayAccess( source_data->getArrayData() )
00199       , source_data->getBox().lower()
00200       , source_data->getBox().upper()
00201       , patch_geom->getXLower()
00202       , patch_geom->getXUpper()
00203       , patch_geom->getDx()
00204       );
00205     printObjectX(plog,source_data->getArrayData());
00206     /* Set diffusion coefficients. */
00207     MultiDimArrayAccess<double,2> t2
00208       = arrayData2ArrayAccess( diffcoef_data->getArrayData(0) );
00209     setArrayDataToConstant(
00210         t2
00211       , diffcoef_data->getGhostBox().lower()
00212       , diffcoef_data->getGhostBox().upper()
00213       , patch_geom->getXLower()
00214       , patch_geom->getXUpper()
00215       , patch_geom->getDx()
00216       , 1.0
00217       );
00218     MultiDimArrayAccess<double,2> t3
00219       = arrayData2ArrayAccess( diffcoef_data->getArrayData(1) );
00220     setArrayDataToConstant(
00221         t3
00222       , diffcoef_data->getGhostBox().lower()
00223       , diffcoef_data->getGhostBox().upper()
00224       , patch_geom->getXLower()
00225       , patch_geom->getXUpper()
00226       , patch_geom->getDx()
00227       , 1.0
00228       );
00229   }
00230   /* Set vector weight. */
00231   setVectorWeightsX( hierarchy, d_weight_persistent );
00232 }
00233 
00234 
00235 
00236 
00237 int JacobiPoissonX::jacobiIteration(
00238   const tbox_Pointer< SAMRAI::hier_PatchHierarchyX > hierarchy
00239 , double relax_coef
00240 ) {
00241 
00242   /*
00243     Jacobi iteration requires temporarily copying solution into scratch space.
00244     In this loop, we allocate the scratch space.  It will be deallocated
00245     before leaving in this function.
00246   */
00247   for ( int ln=0; ln < hierarchy->getNumberLevels(); ++ln ) {
00248     tbox_Pointer<hier_PatchLevelX> level = hierarchy->getPatchLevel(ln);
00249     level->allocatePatchData(d_scalar_scratch);
00250   }
00251 
00252 
00253 
00254   /*
00255     Refine solution into scratch space to fill the fine level
00256     ghost cells in the scratch space.  The interiors will be
00257     filled also, but then they will be overwritten by copying
00258     the actual fine grid solution into the fine coverage area.
00259 
00260     1. Create a refine algorithm.
00261 
00262     2. Associate a refinement operator with source-desination pairs of data
00263     (specified by indices).
00264 
00265     3. Cycle through the levels and refine the data.
00266   */
00267   // 1.
00268   xfer_RefineAlgorithmX refiner;
00269   // 2.
00270   {
00271 
00272     tbox_Pointer<geom_CartesianCellDoubleConservativeLinearRefineX>
00273       conservative_linear_refine_operator
00274       = new geom_CartesianCellDoubleConservativeLinearRefineX;
00275     assert( conservative_linear_refine_operator->getOperatorName() == "CONSERVATIVE_LINEAR_REFINE" );
00276 
00277     tbox_Pointer<geom_CartesianCellDoubleLinearRefineX>
00278       linear_refine_operator
00279       = new geom_CartesianCellDoubleLinearRefineX;
00280     assert( linear_refine_operator->getOperatorName() == "LINEAR_REFINE" );
00281 
00282     refiner.registerRefine( d_scalar_scratch, d_scalar_persistent
00283                           , d_scalar_scratch, conservative_linear_refine_operator );
00284   }
00285   // 3.
00286   for ( int src_ln=0; src_ln < hierarchy->getNumberLevels()-1; ++src_ln ) {
00287     int dst_ln = src_ln+1;
00288     /*
00289       Source level is the the current level (number ln).
00290       Destination level is the next finer level (ln+1).
00291     */
00292     tbox_Pointer<hier_PatchLevelX> src_level
00293       = hierarchy->getPatchLevel(src_ln);
00294     tbox_Pointer<hier_PatchLevelX> dst_level
00295       = hierarchy->getPatchLevel(dst_ln);
00296     assert(src_level);
00297     assert(dst_level);
00298 
00299     {
00300       // refiner.createSchedule( dst_level , src_level )->fillData(0.0);
00301       // tbox_Pointer<xfer_RefineScheduleX> refine_schedule =
00302       // refiner.createSchedule( dst_level , src_level );
00303       tbox_Pointer<xfer_RefineScheduleX> refine_schedule =
00304         refiner.createSchedule( dst_level, src_ln, hierarchy );
00305         // refiner.createSchedule( dst_level, src_level, hierarchy );
00306       refine_schedule->fillData(0.0);
00307     }
00308 
00309   }
00310 
00311 
00312 
00313   /*
00314     Copy solution so that the interiors of scratch patches get set
00315     to the current solution.
00316   */
00317   for ( int ln=0; ln < hierarchy->getNumberLevels(); ++ln ) {
00318     tbox_Pointer<hier_PatchLevelX> level = hierarchy->getPatchLevel(ln);
00319 
00320     for ( hier_PatchLevelX::Iterator p(level); p; p++ ) {
00321       tbox_Pointer<hier_PatchX> patch = level->getPatch(p());
00322 
00323       /* Get the patch data. */
00324       tbox_Pointer<hier_PatchDataX> current_solution
00325         = patch->getPatchData(d_scalar_persistent);
00326       tbox_Pointer<hier_PatchDataX> scratch_solution
00327         = patch->getPatchData(d_scalar_scratch);
00328 
00329       /* Copy the current solution into the scratch space. */
00330       // plog << "current before copying: level " << ln << " patch " << p() ;
00331       // printObjectX( plog, (hier_PatchDataX&)current_solution );
00332       // printObjectX( plog, current_solution->getArrayData() );
00333       scratch_solution->copy(*current_solution);
00334       // plog << "scratch after copying: level " << ln << " patch " << p() ;
00335       // printObjectX( plog, (hier_PatchDataX&)scratch_solution );
00336       // printObjectX( plog, scratch_solution->getArrayData() );
00337 
00338     }
00339 
00340   }
00341 
00342 
00343 
00344 
00345   /*
00346     Perform jacobi relaxation on all levels, all patches.
00347   */
00348   for ( int ln = hierarchy->getFinestLevelNumber(); ln >= 0; ln -- ) {
00349 
00350     tbox_Pointer<hier_PatchLevelX> level = hierarchy->getPatchLevel(ln);
00351 
00352     for ( hier_PatchLevelX::Iterator p(level); p; p++ ) {
00353 
00354       tbox_Pointer<hier_PatchX> patch = level->getPatch(p());
00355 
00356       /*
00357         Get information about the patch geometry.
00358       */
00359       tbox_Pointer<geom_CartesianPatchGeometryX> patch_geom
00360         = patch->getPatchGeometry();
00361 
00362       /*
00363         Get the patch data.
00364       */
00365       tbox_Pointer<pdat_CellDataX<double> > current_solution
00366         = patch->getPatchData(d_scalar_persistent);
00367       tbox_Pointer<pdat_CellDataX<double> > scratch_solution
00368         = patch->getPatchData(d_scalar_scratch);
00369       tbox_Pointer<pdat_CellDataX<double> > source
00370         = patch->getPatchData(d_source_persistent);
00371 
00372       if ( ln == 0 ) {
00373         /*
00374           Set physical boundary condition in the scratch space
00375           for coarsest level.
00376         */
00377         jacobiPoissonBoundary(
00378             arrayData2ArrayAccess( scratch_solution->getArrayData(), 0 )
00379           , current_solution->getBox().lower()
00380           , current_solution->getBox().upper()
00381           , patch_geom->getXLower()
00382           , patch_geom->getXUpper()
00383           , patch_geom->getDx()
00384           );
00385       }
00386       if( 0 ) jacobiPoissonExactBC(
00387           arrayData2ArrayAccess( scratch_solution->getArrayData() )
00388         , current_solution->getBox().lower()
00389         , current_solution->getBox().upper()
00390         , patch_geom->getXLower()
00391         , patch_geom->getXUpper()
00392         , patch_geom->getDx()
00393         );
00394 
00395       /*
00396         Relax solution in scratch and assign result to the solution.
00397       */
00398       // plog << "scratch before: level " << ln << " patch " << p() ;
00399       // printObjectX( plog, (hier_PatchDataX&)scratch_solution );
00400       // printObjectX( plog, scratch_solution->getArrayData() );
00401       // plog << "current before: level " << ln << " patch " << p() ;
00402       // printObjectX( plog, (hier_PatchDataX&)current_solution );
00403       // printObjectX( plog, current_solution->getArrayData() );
00404       jacobiPoissonSmooth(
00405           relax_coef
00406         , arrayData2ArrayAccess( current_solution->getArrayData(), 0 )
00407         , arrayData2ArrayAccess( scratch_solution->getArrayData(), 0 )
00408         , arrayData2ArrayAccess( source->getArrayData(), 0 )
00409         , current_solution->getBox().lower()
00410         , current_solution->getBox().upper()
00411         , patch_geom->getXLower()
00412         , patch_geom->getXUpper()
00413         , patch_geom->getDx()
00414         );
00415       // plog << "current after: level " << ln << " patch " << p() ;
00416       // printObjectX( plog, (hier_PatchDataX&)current_solution );
00417       // printObjectX( plog, current_solution->getArrayData() );
00418 
00419     }
00420 
00421   }
00422 
00423 
00424 
00425   /*
00426     The computational stencils will be extending into the ghost cells
00427     of the scratch space, so for each level (except the finest) get
00428     the data into the ghost cells using data from the next finest level.
00429 
00430     1. Create a coarsen algorithm to coarsen data on fine patch interiors.
00431 
00432     2.
00433     Associate a coarsening operator with source-desination pairs of data
00434     (specified by indices).  The indices are the same because we transfer
00435     data from between levels, but it is the same variable in the same
00436     context.
00437 
00438     3. Cycle through the levels and coarsen the data.
00439   */
00440 
00441   // 1.
00442   xfer_CoarsenAlgorithmX coarsener;
00443   // 2.
00444   {
00445     tbox_Pointer<geom_CartesianCellDoubleWeightedAverageX>
00446       conservative_coarsen_operator
00447       = new geom_CartesianCellDoubleWeightedAverageX;
00448     coarsener.registerCoarsen(
00449                                d_scalar_persistent
00450                              , d_scalar_persistent
00451                              , conservative_coarsen_operator
00452                              );
00453   }
00454   // 3.
00455   for ( int ln = hierarchy->getFinestLevelNumber(); ln > 0; --ln ) {
00456     if(1) coarsener.createSchedule(
00457         hierarchy->getPatchLevel(ln-1)
00458       , hierarchy->getPatchLevel(ln)
00459       )->coarsenData();
00460   }
00461 
00462   /*
00463     Deallocate scratch data.
00464   */
00465   for ( int ln=0; ln < hierarchy->getNumberLevels(); ++ln ) {
00466     tbox_Pointer<hier_PatchLevelX> level = hierarchy->getPatchLevel(ln);
00467     level->deallocatePatchData(d_scalar_scratch);
00468   }
00469 
00470   return 0;
00471 }
00472 
00473 
00474 
00475 
00476 
00477 
00478 int JacobiPoissonX::registerVariablesWithPlotter(
00479   tbox_Pointer<plot_CartesianVizamraiDataWriterX> viz_writer
00480 ) const {
00481 
00482   /*
00483     Register variables with plotter.
00484   */
00485   viz_writer->registerPlotVariable("Computed solution", d_scalar_persistent);
00486   viz_writer->registerPlotVariable("Exact solution", d_exact_persistent);
00487   viz_writer->registerDerivedPlotVariable("Error"
00488     , (plot_VizamraiDerivedDataStrategyX *)this );
00489   viz_writer->registerPlotVariable("Poisson source", d_source_persistent);
00490   viz_writer->registerPlotVariable("Vector weight", d_weight_persistent);
00491 
00492   return 0;
00493 }
00494 
00495 
00496 
00497 
00498 void JacobiPoissonX::writeDerivedDataToStream(
00499   tbox_AbstractStream &stream
00500 , const hier_PatchX &patch
00501 , int patch_level_number
00502 , const hier_BoxX &region
00503 , const string &variable_name
00504 , int plot_type
00505 ) {
00506 
00507   assert(plot_type==1); // Writing float data.
00508   assert( variable_name == "Error" );
00509   /*
00510     Get the patch data.
00511   */
00512   tbox_Pointer<pdat_CellDataX<double> > current_solution
00513     = patch.getPatchData(d_scalar_persistent);
00514   tbox_Pointer<pdat_CellDataX<double> > exact_solution
00515     = patch.getPatchData(d_exact_persistent);
00516 
00517   {
00518     ConstMultiDimArrayAccess<double,2> ex =
00519       arrayData2ArrayAccess( exact_solution->getArrayData() );
00520     ConstMultiDimArrayAccess<double,2> co =
00521       arrayData2ArrayAccess( current_solution->getArrayData() );
00522     const int *lower = region.lower();
00523     const int *upper = region.upper();
00524     for ( int j=lower[1]; j<=upper[1]; ++j ) {
00525       for ( int i=lower[0]; i<=upper[0]; ++i ) {
00526         float diff = ( co(j,i) - ex(j,i) );
00527         stream.pack(&diff);
00528       }
00529     }
00530   }
00531   return;
00532 }
00533 
00534 
00535 
00536 
00537 
00538 int JacobiPoissonX::computeError(
00539   const tbox_Pointer< SAMRAI::hier_PatchHierarchyX > hierarchy
00540 , double *l2norm
00541 , double *linorm
00542 ) const {
00543 
00544   /*
00545     Compute error on all levels, all patches.
00546     This can be done using the data operator classes in SAMRAI.
00547     But we do it manually here to combine the operations into
00548     the least number of loopings and to illustrate how it is done.
00549   */
00550   double diff=0, difflinf=0, diffl2=0, weightsum=0;
00551   for ( int ln = hierarchy->getFinestLevelNumber(); ln >= 0; ln -- ) {
00552     tbox_Pointer<hier_PatchLevelX> level = hierarchy->getPatchLevel(ln);
00553 
00554     for ( hier_PatchLevelX::Iterator p(level); p; p++ ) {
00555       tbox_Pointer<hier_PatchX> patch = level->getPatch(p());
00556 
00557       /*
00558         Get the patch data.
00559       */
00560       tbox_Pointer<pdat_CellDataX<double> > current_solution
00561         = patch->getPatchData(d_scalar_persistent);
00562       tbox_Pointer<pdat_CellDataX<double> > exact_solution
00563         = patch->getPatchData(d_exact_persistent);
00564       tbox_Pointer<pdat_CellDataX<double> > weight
00565         = patch->getPatchData(d_weight_persistent);
00566 
00567       {
00568         ConstMultiDimArrayAccess<double,2> ex =
00569           arrayData2ArrayAccess( exact_solution->getArrayData() );
00570         ConstMultiDimArrayAccess<double,2> co =
00571           arrayData2ArrayAccess( current_solution->getArrayData() );
00572         ConstMultiDimArrayAccess<double,2> wt =
00573           arrayData2ArrayAccess( weight->getArrayData() );
00574         const int *lower = current_solution->getBox().lower();
00575         const int *upper = current_solution->getBox().upper();
00576         for ( int j=lower[1]; j<=upper[1]; ++j ) {
00577           for ( int i=lower[0]; i<=upper[0]; ++i ) {
00578             /*
00579               Disregard zero weights in error computations
00580               because they are on coarse grids covered by finer grids.
00581             */
00582             if ( wt(j,i) > 0 ) {
00583               diff = fabs( co(j,i) - ex(j,i) );
00584               if ( difflinf < diff ) difflinf = diff;
00585               diffl2 += wt(j,i)*diff*diff;
00586               weightsum += wt(j,i);
00587             }
00588           }
00589         }
00590       }
00591 
00592     }
00593 
00594   }
00595   diffl2 = sqrt( diffl2/weightsum );
00596 
00597   *l2norm = diffl2;
00598   *linorm = difflinf;
00599 
00600   return 0;
00601 }

Generated on Wed Apr 17 12:51:44 2002 for samtut by doxygen1.2.12 written by Dimitri van Heesch, © 1997-2001