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); * * } */ }
/// <summary> /// Solution of the system /// <see cref="LaplaceMtx"/>*<see cref="T"/> + <see cref="LaplaceAffine"/> = <see cref="RHS"/> /// using the modular solver framework. /// </summary> private void ExperimentalSolve(out double mintime, out double maxtime, out bool Converged, out int NoOfIter) { using (var tr = new FuncTrace()) { int p = this.T.Basis.Degree; var MgSeq = this.MultigridSequence; mintime = double.MaxValue; maxtime = 0; Converged = false; NoOfIter = int.MaxValue; Console.WriteLine("Construction of Multigrid basis..."); Stopwatch mgBasis = new Stopwatch(); mgBasis.Start(); AggregationGridBasis[][] AggBasis; using (new BlockTrace("Aggregation_basis_init", tr)) { AggBasis = AggregationGridBasis.CreateSequence(MgSeq, new Basis[] { this.T.Basis }); } mgBasis.Stop(); Console.WriteLine("done. (" + mgBasis.Elapsed.TotalSeconds + " sec)"); //foreach (int sz in new int[] { 1000, 2000, 5000, 10000, 20000 }) { // base.Control.TargetBlockSize = sz; for (int irun = 0; irun < base.Control.NoOfSolverRuns; irun++) { Stopwatch stw = new Stopwatch(); stw.Reset(); stw.Start(); Console.WriteLine("Setting up multigrid operator..."); var mgsetup = new Stopwatch(); mgsetup.Start(); var MultigridOp = new MultigridOperator(AggBasis, this.T.Mapping, this.LaplaceMtx, null, MgConfig); mgsetup.Stop(); Console.WriteLine("done. (" + mgsetup.Elapsed.TotalSeconds + " sec)"); Console.WriteLine("Setting up solver..."); var solverSetup = new Stopwatch(); solverSetup.Start(); ISolverSmootherTemplate solver; switch (base.Control.solver_name) { case SolverCodes.exp_direct: solver = new DirectSolver() { WhichSolver = DirectSolver._whichSolver.PARDISO }; break; case SolverCodes.exp_direct_lapack: solver = new DirectSolver() { WhichSolver = DirectSolver._whichSolver.Lapack }; break; case SolverCodes.exp_softpcg_schwarz_directcoarse: { double LL = this.LaplaceMtx._RowPartitioning.LocalLength; int NoOfBlocks = (int)Math.Max(1, Math.Round(LL / (double)this.Control.TargetBlockSize)); Console.WriteLine("Additive Schwarz w. direct coarse, No of blocks: " + NoOfBlocks.MPISum()); solver = new SoftPCG() { m_MaxIterations = 50000, m_Tolerance = 1.0e-10, Precond = new Schwarz() { m_MaxIterations = 1, //CoarseSolver = new GenericRestriction() { // CoarserLevelSolver = new GenericRestriction() { CoarseSolver = new DirectSolver() { WhichSolver = DirectSolver._whichSolver.PARDISO // } //} }, m_BlockingStrategy = new Schwarz.METISBlockingStrategy() { NoOfPartsPerProcess = NoOfBlocks }, Overlap = 1, } }; break; } case SolverCodes.exp_softpcg_schwarz: { double LL = this.LaplaceMtx._RowPartitioning.LocalLength; int NoOfBlocks = (int)Math.Max(1, Math.Round(LL / (double)this.Control.TargetBlockSize)); Console.WriteLine("Additive Schwarz, No of blocks: " + NoOfBlocks.MPISum()); solver = new SoftPCG() { m_MaxIterations = 50000, m_Tolerance = 1.0e-10, Precond = new Schwarz() { m_MaxIterations = 1, CoarseSolver = null, m_BlockingStrategy = new Schwarz.METISBlockingStrategy { NoOfPartsPerProcess = NoOfBlocks }, Overlap = 1 } }; break; } case SolverCodes.exp_softpcg_mg: solver = MultilevelSchwarz(MultigridOp); break; case SolverCodes.exp_Kcycle_schwarz: solver = KcycleMultiSchwarz(MultigridOp); break; default: throw new ApplicationException("unknown solver: " + this.Control.solver_name); } T.Clear(); T.AccLaidBack(1.0, Tex); ConvergenceObserver CO = null; //CO = new ConvergenceObserver(MultigridOp, null, T.CoordinateVector.ToArray()); //CO.TecplotOut = "oasch"; if (solver is ISolverWithCallback) { if (CO == null) { ((ISolverWithCallback)solver).IterationCallback = delegate (int iter, double[] xI, double[] rI, MultigridOperator mgOp) { double l2_RES = rI.L2NormPow2().MPISum().Sqrt(); double[] xRef = new double[xI.Length]; MultigridOp.TransformSolInto(T.CoordinateVector, xRef); double l2_ERR = GenericBlas.L2DistPow2(xI, xRef).MPISum().Sqrt(); Console.WriteLine("Iter: {0}\tRes: {1:0.##E-00}\tErr: {2:0.##E-00}\tRunt: {3:0.##E-00}", iter, l2_RES, l2_ERR, stw.Elapsed.TotalSeconds); //Tjac.CoordinatesAsVector.SetV(xI); //Residual.CoordinatesAsVector.SetV(rI); //PlotCurrentState(iter, new TimestepNumber(iter), 3); }; } else { ((ISolverWithCallback)solver).IterationCallback = CO.IterationCallback; } } using (new BlockTrace("Solver_Init", tr)) { solver.Init(MultigridOp); } solverSetup.Stop(); Console.WriteLine("done. (" + solverSetup.Elapsed.TotalSeconds + " sec)"); Console.WriteLine("Running solver..."); var solverIteration = new Stopwatch(); solverIteration.Start(); double[] T2 = this.T.CoordinateVector.ToArray(); using (new BlockTrace("Solver_Run", tr)) { solver.ResetStat(); T2.Clear(); var RHSvec = RHS.CoordinateVector.ToArray(); BLAS.daxpy(RHSvec.Length, -1.0, this.LaplaceAffine, 1, RHSvec, 1); MultigridOp.UseSolver(solver, T2, RHSvec); T.CoordinateVector.SetV(T2); } solverIteration.Stop(); Console.WriteLine("done. (" + solverIteration.Elapsed.TotalSeconds + " sec)"); Console.WriteLine("Pardiso phase 11: " + ilPSP.LinSolvers.PARDISO.PARDISOSolver.Phase_11.Elapsed.TotalSeconds); Console.WriteLine("Pardiso phase 22: " + ilPSP.LinSolvers.PARDISO.PARDISOSolver.Phase_22.Elapsed.TotalSeconds); Console.WriteLine("Pardiso phase 33: " + ilPSP.LinSolvers.PARDISO.PARDISOSolver.Phase_33.Elapsed.TotalSeconds); // time measurement, statistics stw.Stop(); mintime = Math.Min(stw.Elapsed.TotalSeconds, mintime); maxtime = Math.Max(stw.Elapsed.TotalSeconds, maxtime); Converged = solver.Converged; NoOfIter = solver.ThisLevelIterations; if (CO != null) CO.PlotTrend(true, true, true); } } }