private void ConsistencyTest() { // consistency test on the original matrix // ----------------------------------------------------- this.residual.Clear(); double[] RHSvec = this.GetRHS(); this.residual.CoordinateVector.SetV(RHSvec, 1.0); this.Op_Matrix.SpMV(-1.0, this.u.CoordinateVector, 1.0, residual.CoordinateVector); double residual_L2Norm = this.residual.L2Norm(); Console.WriteLine("Residual norm: " + residual_L2Norm); Assert.LessOrEqual(residual_L2Norm, 1.0e-8); // consistency test on the multigrid // -------------------------------------------------------------------- AggregationGridBasis[][] XAggB = AggregationGridBasis.CreateSequence(base.MultigridSequence, u.Mapping.BasisS); XAggB.UpdateXdgAggregationBasis(this.Op_Agglomeration); int p = this.u.Basis.Degree; var MultigridOp = new MultigridOperator(XAggB, this.u.Mapping, this.Op_Matrix, this.Op_mass.GetMassMatrix(new UnsetteledCoordinateMapping(this.u.Basis), false), new MultigridOperator.ChangeOfBasisConfig[][] { new MultigridOperator.ChangeOfBasisConfig[] { new MultigridOperator.ChangeOfBasisConfig() { VarIndex = new int[] { 0 }, mode = MultigridOperator.Mode.Eye, Degree = u.Basis.Degree } } }); double[] mgSolVec = new double[MultigridOp.Mapping.LocalLength]; double[] mgRhsVec = new double[MultigridOp.Mapping.LocalLength]; MultigridOp.TransformSolInto(this.u.CoordinateVector, mgSolVec); MultigridOp.TransformRhsInto(RHSvec, mgRhsVec); MgConsistencyTestRec(MultigridOp, mgSolVec, mgRhsVec); // /* * { * int Jagg1 = MgSeq[1].NoOfAggregateCells; * MultigridOperator MgOp0 = MultigridOp; * MultigridOperator MgOp1 = MultigridOp.CoarserLevel; * MultigridMapping Map0 = MgOp0.Mapping; * MultigridMapping Map1 = MgOp1.Mapping; * * * double[] V0 = new double[MgOp0.Mapping.LocalLength]; * double[] V1 = new double[MgOp1.Mapping.LocalLength]; * * for(int j = 0; j < Jagg1; j++) { * int idx = Map1.LocalUniqueIndex(0, j, 0); * V1[idx] = j; * } * * MgOp1.Prolongate(1.0, V0, 0.0, V1); * * XDGField Marker = new XDGField(this.u.Basis, "Tracker"); * * MgOp0.TransformSolFrom(Marker.CoordinatesAsVector, V0); * this.Op_Agglomeration.Extrapolate(Marker.CoordinatesAsVector, Marker.Mapping); * * Tecplot.PlotFields(new DGField[] { Marker }, "Tracker", "Tracker", 0.0, 5); * * } */ }
private void ExperimentalSolver(out double mintime, out double maxtime, out bool Converged, out int NoOfIter, out int DOFs) { using (var tr = new FuncTrace()) { mintime = double.MaxValue; maxtime = 0; Converged = false; NoOfIter = int.MaxValue; DOFs = 0; AggregationGridBasis[][] XAggB; using (new BlockTrace("Aggregation_basis_init", tr)) { XAggB = AggregationGridBasis.CreateSequence(base.MultigridSequence, u.Mapping.BasisS); } XAggB.UpdateXdgAggregationBasis(this.Op_Agglomeration); var MassMatrix = this.Op_mass.GetMassMatrix(this.u.Mapping, new double[] { 1.0 }, false, this.LsTrk.SpeciesIdS.ToArray()); double[] _RHSvec = this.GetRHS(); Stopwatch stw = new Stopwatch(); stw.Reset(); stw.Start(); Console.WriteLine("Setting up multigrid operator..."); int p = this.u.Basis.Degree; var MultigridOp = new MultigridOperator(XAggB, this.u.Mapping, this.Op_Matrix, this.Op_mass.GetMassMatrix(new UnsetteledCoordinateMapping(this.u.Basis), false), OpConfig); Assert.True(MultigridOp != null); int L = MultigridOp.Mapping.LocalLength; DOFs = MultigridOp.Mapping.TotalLength; double[] RHSvec = new double[L]; MultigridOp.TransformRhsInto(_RHSvec, RHSvec); ISolverSmootherTemplate exsolver; SolverFactory SF = new SolverFactory(this.Control.NonLinearSolver, this.Control.LinearSolver); List <Action <int, double[], double[], MultigridOperator> > Callbacks = new List <Action <int, double[], double[], MultigridOperator> >(); Callbacks.Add(CustomItCallback); SF.GenerateLinear(out exsolver, MultigridSequence, OpConfig, Callbacks); using (new BlockTrace("Solver_Init", tr)) { exsolver.Init(MultigridOp); } /* * string filename = "XdgPoisson" + this.Grid.SpatialDimension + "p" + this.u.Basis.Degree + "R" + this.Grid.CellPartitioning.TotalLength; * MultigridOp.OperatorMatrix.SaveToTextFileSparse(filename + ".txt"); * RHSvec.SaveToTextFile(filename + "_rhs.txt"); * * var uEx = this.u.CloneAs(); * Op_Agglomeration.ClearAgglomerated(uEx.Mapping); * var CO = new ConvergenceObserver(MultigridOp, MassMatrix, uEx.CoordinateVector.ToArray()); * uEx = null; * CO.TecplotOut = "PoissonConvergence"; * //CO.PlotDecomposition(this.u.CoordinateVector.ToArray()); * * if (exsolver is ISolverWithCallback) { * ((ISolverWithCallback)exsolver).IterationCallback = CO.IterationCallback; * } * //*/ XDGField u2 = u.CloneAs(); using (new BlockTrace("Solver_Run", tr)) { // use solver (on XDG-field 'u2'). u2.Clear(); MultigridOp.UseSolver(exsolver, u2.CoordinateVector, _RHSvec); Console.WriteLine("Solver: {0}, converged? {1}, {2} iterations.", exsolver.GetType().Name, exsolver.Converged, exsolver.ThisLevelIterations); this.Op_Agglomeration.Extrapolate(u2.Mapping); Assert.IsTrue(exsolver.Converged, "Iterative solver did not converge."); } stw.Stop(); mintime = Math.Min(stw.Elapsed.TotalSeconds, mintime); maxtime = Math.Max(stw.Elapsed.TotalSeconds, maxtime); Converged = exsolver.Converged; NoOfIter = exsolver.ThisLevelIterations; // compute error between reference solution and multigrid solver XDGField ErrField = u2.CloneAs(); ErrField.Acc(-1.0, u); double ERR = ErrField.L2Norm(); double RelERR = ERR / u.L2Norm(); Assert.LessOrEqual(RelERR, 1.0e-6, "Result from iterative solver above threshold."); } }