public XDGTestSetup( int p, double AggregationThreshold, int TrackerWidth, MultigridOperator.Mode mumo, XQuadFactoryHelper.MomentFittingVariants momentFittingVariant, ScalarFunction LevSetFunc = null) { // Level set, tracker and XDG basis // ================================ if (LevSetFunc == null) { LevSetFunc = ((_2D)((x, y) => 0.8 * 0.8 - x * x - y * y)).Vectorize(); } LevSet = new LevelSet(new Basis(grid, 2), "LevelSet"); LevSet.Clear(); LevSet.ProjectField(LevSetFunc); LsTrk = new LevelSetTracker(grid, XQuadFactoryHelper.MomentFittingVariants.Classic, TrackerWidth, new string[] { "A", "B" }, LevSet); LsTrk.UpdateTracker(); XB = new XDGBasis(LsTrk, p); XSpatialOperator Dummy = new XSpatialOperator(1, 0, 1, QuadOrderFunc.SumOfMaxDegrees(RoundUp: true), "C1", "u"); //Dummy.EquationComponents["c1"].Add(new Dummy.Commit(); //Tecplot.PlotFields(new DGField[] { LevSet }, "agglo", 0.0, 3); // operator // ======== Debug.Assert(p <= 4); XDGBasis opXB = new XDGBasis(LsTrk, 4); // we want to have a very precise quad rule var map = new UnsetteledCoordinateMapping(opXB); int quadOrder = Dummy.QuadOrderFunction(map.BasisS.Select(bs => bs.Degree).ToArray(), new int[0], map.BasisS.Select(bs => bs.Degree).ToArray()); //agg = new MultiphaseCellAgglomerator(new CutCellMetrics(momentFittingVariant, quadOrder, LsTrk, LsTrk.SpeciesIdS.ToArray()), AggregationThreshold, false); agg = LsTrk.GetAgglomerator(LsTrk.SpeciesIdS.ToArray(), quadOrder, __AgglomerationTreshold: AggregationThreshold); foreach (var S in LsTrk.SpeciesIdS) { Console.WriteLine("Species {0}, no. of agglomerated cells {1} ", LsTrk.GetSpeciesName(S), agg.GetAgglomerator(S).AggInfo.SourceCells.Count()); } // mass matrix factory // =================== // Basis maxB = map.BasisS.ElementAtMax(b => b.Degree); //MassFact = new MassMatrixFactory(maxB, agg); MassFact = LsTrk.GetXDGSpaceMetrics(LsTrk.SpeciesIdS.ToArray(), quadOrder, 1).MassMatrixFactory; // Test field // ========== // set the test field: this is a polynomial function, // but different for each species; On this field, restriction followed by prolongation should be the identity this.Xdg_uTest = new XDGField(this.XB, "uTest"); Dictionary <SpeciesId, double> dumia = new Dictionary <SpeciesId, double>(); int i = 2; foreach (var Spc in LsTrk.SpeciesIdS) { dumia.Add(Spc, i); i -= 1; } SetTestValue(Xdg_uTest, dumia); // dummy operator matrix which fits polynomial degree p // ==================================================== Xdg_opMtx = new BlockMsrMatrix(Xdg_uTest.Mapping, Xdg_uTest.Mapping); Xdg_opMtx.AccEyeSp(120.0); // XDG Aggregation BasiseS // ======================= //XAggB = MgSeq.Select(agGrd => new XdgAggregationBasis[] { new XdgAggregationBasis(uTest.Basis, agGrd) }).ToArray(); XAggB = new XdgAggregationBasis[MgSeq.Length][]; var _XAggB = AggregationGridBasis.CreateSequence(MgSeq, Xdg_uTest.Mapping.BasisS); for (int iLevel = 0; iLevel < XAggB.Length; iLevel++) { XAggB[iLevel] = new[] { (XdgAggregationBasis)(_XAggB[iLevel][0]) }; XAggB[iLevel][0].Update(agg); } // Multigrid Operator // ================== Xdg_opMtx = new BlockMsrMatrix(Xdg_uTest.Mapping, Xdg_uTest.Mapping); Xdg_opMtx.AccEyeSp(120.0); XdgMultigridOp = new MultigridOperator(XAggB, Xdg_uTest.Mapping, Xdg_opMtx, MassFact.GetMassMatrix(Xdg_uTest.Mapping, false), new MultigridOperator.ChangeOfBasisConfig[][] { new MultigridOperator.ChangeOfBasisConfig[] { new MultigridOperator.ChangeOfBasisConfig() { VarIndex = new int[] { 0 }, mode = mumo, Degree = p } } }); }
public static void XDG_ProlongationTest( [Values(0, 1, 2, 3)] int p, [Values(0.0, 0.3)] double AggregationThreshold, [Values(0, 1)] int TrackerWidth, [Values(MultigridOperator.Mode.Eye, MultigridOperator.Mode.IdMass)] MultigridOperator.Mode mode) { XQuadFactoryHelper.MomentFittingVariants variant = XQuadFactoryHelper.MomentFittingVariants.OneStepGaussAndStokes; var xt = new XDGTestSetup(p, AggregationThreshold, TrackerWidth, MultigridOperator.Mode.Eye, variant //, ((Func<double[], double>)(X => X[0] + 0.75)).Vectorize() ); int Jup = grid.Cells.NoOfLocalUpdatedCells; Random rnd = new Random(); int Ltop = xt.XdgMultigridOp.Mapping.LocalLength; // Number of DOF's on top multigrid level. double[] RndVec = Ltop.ForLoop(i => rnd.NextDouble()); double[] NoJmpVec = new double[Ltop]; for (int iLevel = 0; iLevel < MgSeq.Length - 1; iLevel++) { XDG_Recursive(0, iLevel, xt.XdgMultigridOp, RndVec, NoJmpVec); // restrict RndVec downt to level 'iLevel', and back up // right now, the XDG field defined by 'NoJmpVec' should be a member // of the aggregated XDG space on level 'iLevel'; // so, there should be no inter-element jumps on the fine level, for each aggregated cell. // Let's test that! XDGField Test = new XDGField(xt.XB, "Test"); xt.XdgMultigridOp.TransformSolFrom(Test.CoordinateVector, NoJmpVec); //xt.agg.Extrapolate(Test.Mapping); var aggGrd = MgSeq[iLevel]; foreach (var spc in xt.LsTrk.SpeciesIdS) { var Test_spc = Test.GetSpeciesShadowField(spc); var SpcMask = xt.LsTrk.Regions.GetSpeciesMask(spc); BitArray AggSourceBitmask = xt.agg.GetAgglomerator(spc).AggInfo.SourceCells.GetBitMask(); double Err = 0; for (int jagg = 0; jagg < aggGrd.iLogicalCells.NoOfLocalUpdatedCells; jagg++) { BitArray CompCellMask = new BitArray(Jup); foreach (int jCell in aggGrd.iLogicalCells.AggregateCellToParts[jagg]) { if (!AggSourceBitmask[jCell]) { CompCellMask[jCell] = true; } } SubGrid CompCellSubGrid = new SubGrid((new CellMask(grid, CompCellMask)).Intersect(SpcMask)); Err += JumpNorm(Test_spc, CompCellSubGrid.InnerEdgesMask).Pow2(); } Console.WriteLine("prolongation jump test (level {0}, species {2}): {1}", iLevel, Err, xt.LsTrk.GetSpeciesName(spc)); Assert.LessOrEqual(Err, 1.0e-8); } } }
public static void XDG_MatrixPolynomialRestAndPrlgTest_2( [Values(0, 1, 2, 3)] int p, [Values(0.0, 0.3)] double AggregationThreshold, [Values(0, 1)] int TrackerWidth, [Values(MultigridOperator.Mode.Eye, MultigridOperator.Mode.IdMass)] MultigridOperator.Mode mode) { if (AggregationThreshold < 0.1 && p >= 3 && mode == MultigridOperator.Mode.IdMass) { // this test combination is not supposed to work: // without agglomeration, for high p, the mass matrix may be indefinite in small cut-cells // => Cholesky decomposition on mass matrix fails, i.e. 'mode == IdMass' cannot succseed. return; } XQuadFactoryHelper.MomentFittingVariants variant = XQuadFactoryHelper.MomentFittingVariants.OneStepGaussAndStokes; var xt = new XDGTestSetup(p, AggregationThreshold, TrackerWidth, mode, variant); // Restriction & prolongation together with orthonormalization // ----------------------------------------------------------- for (var mgop = xt.MultigridOp; mgop != null; mgop = mgop.CoarserLevel) { var Itself = mgop.Mapping.FromOtherLevelMatrix(mgop.Mapping); Itself.AccEyeSp(-1.0); double Itslef_Norm = Itself.InfNorm(); //Console.WriteLine("Level {0}, Restriction onto itself {1}", mgm.LevelIndex, Itslef_Norm); Assert.LessOrEqual(Itslef_Norm, 1.0e-8); } { // test change of basis on top level XDGField uTestRnd = new XDGField(xt.XB); Random rnd = new Random(); for (int i = 0; i < uTestRnd.CoordinateVector.Count; i++) { uTestRnd.CoordinateVector[i] = rnd.NextDouble(); } xt.agg.ClearAgglomerated(uTestRnd.CoordinateVector, uTestRnd.Mapping); // perform change of basis on top level ... int Ltop = xt.MultigridOp.Mapping.LocalLength; double[] uTest_Fine = new double[Ltop]; xt.MultigridOp.TransformSolInto(uTestRnd.CoordinateVector, uTest_Fine); // .. and back XDGField uError2 = uTestRnd.CloneAs(); uError2.Clear(); xt.MultigridOp.TransformSolFrom(uError2.CoordinateVector, uTest_Fine); // compare: uError2.Acc(-1.0, uTestRnd); double NORM_uError = uError2.L2Norm(); // output Console.WriteLine("Top level change of basis error: {0}", NORM_uError); Assert.LessOrEqual(NORM_uError, 1.0e-8); } { // perform change of basis on top level int Ltop = xt.MultigridOp.Mapping.LocalLength; double[] uTest_Fine = new double[Ltop]; xt.MultigridOp.TransformSolInto(xt.uTest.CoordinateVector, uTest_Fine); // check for each level of the multigrid operator... for (int iLevel = 0; iLevel < MgSeq.Count() - 1; iLevel++) { double[] uTest_Prolonged = new double[Ltop]; XDG_Recursive(0, iLevel, xt.MultigridOp, uTest_Fine, uTest_Prolonged); XDGField uError = xt.uTest.CloneAs(); uError.Clear(); xt.MultigridOp.TransformSolFrom(uError.CoordinateVector, uTest_Prolonged); xt.agg.Extrapolate(uError.Mapping); uError.Acc(-1.0, xt.uTest); double NORM_uError = uError.L2Norm(); Console.WriteLine("Rest/Prlg error, level {0}: {1}", iLevel, NORM_uError); Assert.LessOrEqual(NORM_uError, 1.0e-8); } } }