deal.II version GIT relicensing-2173-gae8fc9d14b 2024-11-24 06:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Crystal_Growth_Phase_Field_Model.h
Go to the documentation of this file.
1
117 *  
118 *   #include "PhaseFieldSolver.h"
119 *  
120 *   void InitialValues::vector_value(const Point<2> &p,
121 *   Vector<double> & values) const
122 *   {
123 *   values(0)= 0.0; //Initial p value of domain
124 *   values(1)= 0.2; //Initial temperature of domain
125 *   }
126 * @endcode
127
128
129<a name="ann-PhaseFieldSolver.cpp"></a>
130<h1>Annotated version of PhaseFieldSolver.cpp</h1>
131 *
132 *
133 *
134 *
135 * @code
136 *   /* -----------------------------------------------------------------------------
137 *   *
138 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
139 *   * Copyright (C) 2024 by Umair Hussain
140 *   *
141 *   * This file is part of the deal.II code gallery.
142 *   *
143 *   * -----------------------------------------------------------------------------
144 *   */
145 *  
146 *   #include "PhaseFieldSolver.h"
147 *  
148 *   PhaseFieldSolver::PhaseFieldSolver()
149 *   : mpi_communicator(MPI_COMM_WORLD)
150 *   , n_mpi_processes(Utilities::MPI::n_mpi_processes(mpi_communicator))
151 *   , this_mpi_process(Utilities::MPI::this_mpi_process(mpi_communicator))
152 *   , pcout(std::cout, (this_mpi_process == 0))
153 *   , fe(FE_Q<2>(1), 2)
154 *   , dof_handler(triangulation)
155 *   , time(0.0)
156 *   , final_time(1.)
157 *   , time_step(.0002)
158 *   , theta(0.5)
159 *   , epsilon(0.01)
160 *   , tau(0.0003)
161 *   , gamma(10.)
162 *   , latent_heat(1.4)
163 *   , alpha(0.9)
164 *   , t_eq(1.)
165 *   , a(0.01)
166 *   {}
167 * @endcode
168
169
170<a name="ann-PhaseFieldSolver.h"></a>
171<h1>Annotated version of PhaseFieldSolver.h</h1>
172 *
173 *
174 *
175 *
176
177 *
178 *
179 * @code
180 *   /* -----------------------------------------------------------------------------
181 *   *
182 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
183 *   * Copyright (C) 2024 by Umair Hussain
184 *   *
185 *   * This file is part of the deal.II code gallery.
186 *   *
187 *   * -----------------------------------------------------------------------------
188 *   */
189 *  
190 *   #ifndef KOBAYASHI_PARALLEL_PHASEFIELDSOLVER_H
191 *   #define KOBAYASHI_PARALLEL_PHASEFIELDSOLVER_H
192 *  
193 *   #include <deal.II/base/quadrature_lib.h>
194 *   #include <deal.II/base/function.h>
195 *   #include <deal.II/base/utilities.h>
196 *   #include <deal.II/lac/vector.h>
197 *   #include <deal.II/lac/full_matrix.h>
198 *   #include <deal.II/lac/sparse_matrix.h>
199 *   #include <deal.II/lac/sparse_direct.h>
200 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
201 *   #include <deal.II/lac/solver_cg.h>
202 *   #include <deal.II/lac/precondition.h>
203 *   #include <deal.II/lac/affine_constraints.h>
204 *   #include <deal.II/grid/tria.h>
205 *   #include <deal.II/grid/grid_generator.h>
206 *   #include <deal.II/dofs/dof_handler.h>
207 *   #include <deal.II/dofs/dof_tools.h>
208 *   #include <deal.II/fe/fe_q.h>
209 *   #include <deal.II/fe/fe_values.h>
210 *   #include <deal.II/fe/fe_system.h>
211 *   #include <deal.II/numerics/vector_tools.h>
212 *   #include <deal.II/numerics/matrix_tools.h>
213 *   #include <deal.II/numerics/data_out.h>
214 *   #include <deal.II/grid/grid_in.h>
215 *  
216 * @endcode
217 *
218 * For Parallel Computation
219 *
220 * @code
221 *   #include <deal.II/base/conditional_ostream.h>
222 *   #include <deal.II/base/mpi.h>
223 *   #include <deal.II/lac/petsc_vector.h>
224 *   #include <deal.II/lac/petsc_sparse_matrix.h>
225 *   #include <deal.II/lac/petsc_solver.h>
226 *   #include <deal.II/lac/petsc_precondition.h>
227 *   #include <deal.II/grid/grid_tools.h>
228 *   #include <deal.II/dofs/dof_renumbering.h>
229 *  
230 *   #include <deal.II/distributed/solution_transfer.h>
231 *  
232 *   #include <fstream>
233 *   #include <iostream>
234 *  
235 *   using namespace dealii;
236 *  
237 *   class PhaseFieldSolver {
238 *   public:
239 *   PhaseFieldSolver();
240 *   void run();
241 *  
242 *   private:
243 *   void make_grid_and_dofs();
244 *   void assemble_system();
245 *   void solve();
246 *   void output_results(const unsigned int timestep_number) const;
247 *   double compute_residual();
248 *   void applying_bc();
249 *   float get_random_number();
250 *  
251 *   MPI_Comm mpi_communicator;
252 *   const unsigned int n_mpi_processes;
253 *   const unsigned int this_mpi_process;
254 *   ConditionalOStream pcout;
255 *  
257 *   FESystem<2> fe;
258 *   DoFHandler<2> dof_handler;
259 *   GridIn<2> gridin;
260 *  
261 *   PETScWrappers::MPI::SparseMatrix jacobian_matrix;
262 *  
263 *   double time;
264 *   const double final_time, time_step;
265 *   const double theta;
266 *   const double epsilon, tau, gamma, latent_heat, alpha, t_eq, a; //as given in Ref. [1]
267 *  
268 *   PETScWrappers::MPI::Vector conv_solution; //solution vector at last newton-raphson iteration
269 *   PETScWrappers::MPI::Vector old_solution; //solution vector at last time step
270 *   PETScWrappers::MPI::Vector solution_update; //increment in solution or delta solution
271 *   PETScWrappers::MPI::Vector system_rhs; //to store residual
272 *   Vector<double> conv_solution_np, old_solution_np; //creating non parallel vectors to store data for easy access of old solution values by all processes
273 *  
274 *   };
275 *  
276 * @endcode
277 *
278 * Initial values class
279 *
280 * @code
281 *   class InitialValues : public Function<2>
282 *   {
283 *   public:
284 *   InitialValues(): Function<2>(2)
285 *   {}
286 *   virtual void vector_value(const Point<2> &p,
287 *   Vector<double> & value) const override;
288 *   };
289 *  
290 *  
291 *   #endif //KOBAYASHI_PARALLEL_PHASEFIELDSOLVER_H
292 * @endcode
293
294
295<a name="ann-applying_bc.cpp"></a>
296<h1>Annotated version of applying_bc.cpp</h1>
297 *
298 *
299 *
300 *
301
302 *
303 *
304 * @code
305 *   /* -----------------------------------------------------------------------------
306 *   *
307 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
308 *   * Copyright (C) 2024 by Umair Hussain
309 *   *
310 *   * This file is part of the deal.II code gallery.
311 *   *
312 *   * -----------------------------------------------------------------------------
313 *   */
314 *  
315 *   #include "PhaseFieldSolver.h"
316 *  
317 *   void PhaseFieldSolver::applying_bc(){
318 *   FEValuesExtractors::Scalar phase_parameter(0);
319 *   FEValuesExtractors::Scalar temperature(1);
320 *  
321 *   QGauss<2> quadrature_formula(fe.degree + 1);
322 *   FEValues<2> fe_values(fe,
323 *   quadrature_formula,
325 *  
326 *   ComponentMask p_mask = fe.component_mask(phase_parameter);
327 *   ComponentMask t_mask = fe.component_mask(temperature);
328 *  
329 *   std::map<types::global_dof_index,double> boundary_values;
330 *  
331 * @endcode
332 *
333 * Prescribing p=1 at the left face (this will be maintained in the subsequent iterations when zero BC is applied in the Newton-Raphson iterations)
334 *
335 * @code
337 *   1,
339 *   boundary_values,p_mask);
340 *  
341 * @endcode
342 *
343 * To apply the boundary values only to the solution vector without the Jacobian Matrix and RHS Vector
344 *
345 * @code
346 *   for (auto &boundary_value : boundary_values)
347 *   old_solution(boundary_value.first) = boundary_value.second;
348 *  
349 *   old_solution.compress(VectorOperation::insert);
350 *  
351 *   }
352 * @endcode
353
354
355<a name="ann-assemble_system.cpp"></a>
356<h1>Annotated version of assemble_system.cpp</h1>
357 *
358 *
359 *
360 *
361 * @code
362 *   /* -----------------------------------------------------------------------------
363 *   *
364 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
365 *   * Copyright (C) 2024 by Umair Hussain
366 *   *
367 *   * This file is part of the deal.II code gallery.
368 *   *
369 *   * -----------------------------------------------------------------------------
370 *   */
371 *  
372 *   #include "PhaseFieldSolver.h"
373 *   #include <cmath>
374 *  
375 *   void PhaseFieldSolver::assemble_system() {
376 * @endcode
377 *
378 * Separating each variable as a scalar to easily call the respective shape functions
379 *
380 * @code
381 *   FEValuesExtractors::Scalar phase_parameter(0);
382 *   FEValuesExtractors::Scalar temperature(1);
383 *  
384 *   QGauss<2> quadrature_formula(fe.degree + 1);
385 *   FEValues<2> fe_values(fe,
386 *   quadrature_formula,
388 *  
389 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
390 *   const unsigned int n_q_points = quadrature_formula.size();
391 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
392 *  
393 *   Vector<double> cell_rhs(dofs_per_cell);
394 *  
395 * @endcode
396 *
397 * To copy values and gradients of solution from previous iteration
398 * Old Newton iteration
399 *
400 * @code
401 *   std::vector<Tensor<1, 2>> old_newton_solution_gradients_p(n_q_points);
402 *   std::vector<double> old_newton_solution_values_p(n_q_points);
403 *   std::vector<Tensor<1, 2>> old_newton_solution_gradients_t(n_q_points);
404 *   std::vector<double> old_newton_solution_values_t(n_q_points);
405 * @endcode
406 *
407 * Old time step iteration
408 *
409 * @code
410 *   std::vector<Tensor<1, 2>> old_time_solution_gradients_p(n_q_points);
411 *   std::vector<double> old_time_solution_values_p(n_q_points);
412 *   std::vector<Tensor<1, 2>> old_time_solution_gradients_t(n_q_points);
413 *   std::vector<double> old_time_solution_values_t(n_q_points);
414 *  
415 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
416 *   jacobian_matrix.operator=(0.0);
417 *   system_rhs.operator=(0.0);
418 *  
419 *   for (const auto &cell : dof_handler.active_cell_iterators()){
420 *   if (cell->subdomain_id() == this_mpi_process) {
421 *   cell_matrix = 0;
422 *   cell_rhs = 0;
423 *  
424 *   fe_values.reinit(cell);
425 *  
426 * @endcode
427 *
428 * Copying old solution values
429 *
430 * @code
431 *   fe_values[phase_parameter].get_function_values(conv_solution_np,old_newton_solution_values_p);
432 *   fe_values[phase_parameter].get_function_gradients(conv_solution_np,old_newton_solution_gradients_p);
433 *   fe_values[temperature].get_function_values(conv_solution_np,old_newton_solution_values_t);
434 *   fe_values[temperature].get_function_gradients(conv_solution_np,old_newton_solution_gradients_t);
435 *   fe_values[phase_parameter].get_function_values(old_solution_np,old_time_solution_values_p);
436 *   fe_values[phase_parameter].get_function_gradients(old_solution_np,old_time_solution_gradients_p);
437 *   fe_values[temperature].get_function_values(old_solution_np,old_time_solution_values_t);
438 *   fe_values[temperature].get_function_gradients(old_solution_np,old_time_solution_gradients_t);
439 *  
440 *   for (unsigned int q = 0; q < n_q_points; ++q){
441 *   double khi = get_random_number();
442 * @endcode
443 *
444 * Old solution values
445 *
446 * @code
447 *   double p_on = old_newton_solution_values_p[q]; //old newton solution
448 *   auto grad_p_on = old_newton_solution_gradients_p[q];
449 *   double p_ot = old_time_solution_values_p[q]; //old time step solution
450 *   auto grad_p_ot = old_time_solution_gradients_p[q];
451 *   double t_on = old_newton_solution_values_t[q];
452 *   auto grad_t_on = old_newton_solution_gradients_t[q];
453 *   double t_ot = old_time_solution_values_t[q];
454 *   auto grad_t_ot = old_time_solution_gradients_t[q];
455 *   for (unsigned int i = 0; i < dofs_per_cell; ++i){
456 * @endcode
457 *
458 * Shape Functions
459 *
460 * @code
461 *   double psi_i = fe_values[phase_parameter].value(i,q);
462 *   auto grad_psi_i = fe_values[phase_parameter].gradient(i,q);
463 *   double phi_i = fe_values[temperature].value(i,q);
464 *   auto grad_phi_i = fe_values[temperature].gradient(i,q);
465 *   for (unsigned int j = 0; j < dofs_per_cell; ++j){
466 * @endcode
467 *
468 * Shape Functions
469 *
470 * @code
471 *   double psi_j = fe_values[phase_parameter].value(j,q);
472 *   auto grad_psi_j = fe_values[phase_parameter].gradient(j,q);
473 *   double phi_j = fe_values[temperature].value(j,q);
474 *   auto grad_phi_j = fe_values[temperature].gradient(j,q);
475 *  
476 *   double mp = psi_i*(tau*psi_j);
477 *   double kp = grad_psi_i*(std::pow(epsilon,2)*grad_psi_j);
478 *   double m = (alpha/M_PI)*std::atan(gamma*(t_eq - t_on));
479 *   double t1 = (1-p_on)*(p_on-0.5+m);
480 *   double t2 = -(p_on)*(p_on-0.5+m);
481 *   double t3 = (p_on)*(1-p_on);
482 *   double nl_p = psi_i*((t1+t2+t3)*psi_j);
483 * @endcode
484 *
485 * Adding random noise at the interface
486 *
487 * @code
488 *   nl_p -= a*khi*psi_i*((1.0 - 2*(p_on))*psi_j);
489 *   double f1_p= mp + time_step*theta*kp - time_step*theta*nl_p; // doh f1 by doh p (first Jacobian terms)
490 *  
491 *   double t4 = (p_on)*(1-p_on)*(-(alpha*gamma/(M_PI*(1+std::pow((gamma*(t_eq-t_on)),2)))));
492 *   double nl_t = psi_i*(t4*phi_j);
493 *   double f1_t = -time_step*theta*nl_t; // doh f1 by doh t (second Jacobian terms)
494 *  
495 *   double mpt = phi_i*(latent_heat*psi_j);
496 *   double f2_p = -mpt; // doh f2 by doh p (third Jacobian terms)
497 *  
498 *   double mt = phi_i*(phi_j);
499 *   double kt = grad_phi_i*(grad_phi_j);
500 *   double f2_t = mt + time_step*theta*kt; // doh f2 by doh t (fourth Jacobian terms)
501 *  
502 * @endcode
503 *
504 * Assembling Jacobian matrix
505 *
506 * @code
507 *   cell_matrix(i,j) += (f1_p + f1_t + f2_p + f2_t)*fe_values.JxW(q);
508 *  
509 *   }
510 * @endcode
511 *
512 * Finding f1 and f2 at previous iteration for rhs vector
513 *
514 * @code
515 *   double mp_n = psi_i*(tau*p_on);
516 *   double kp_n = grad_psi_i*(std::pow(epsilon,2)*grad_p_on);
517 *   double m_n = (alpha/M_PI)*std::atan(gamma*(t_eq-t_on));
518 *   double nl_n = psi_i*((p_on)*(1-p_on)*(p_on-0.5+m_n));
519 *   double mp_t = psi_i*(tau*p_ot);
520 *   double kp_t = grad_psi_i*(tau*grad_p_ot);
521 *   double m_t = (alpha/M_PI)*std::atan(gamma*(t_eq-t_ot));
522 *   double nl_t = psi_i*(p_ot)*(1-p_ot)*(p_ot-0.5+m_t);
523 * @endcode
524 *
525 * Adding random noise at the interface
526 *
527 * @code
528 *   nl_n -= psi_i*(a*khi*(p_on)*(1-p_on));
529 *   nl_t -= psi_i*(a*khi*(p_ot)*(1-p_ot));
530 *  
531 *   double f1n = mp_n + time_step*theta*kp_n - time_step*theta*nl_n - mp_t + time_step*(1-theta)*kp_t - time_step*(1-theta)*nl_t; //f1 at last newton iteration
532 *  
533 *   double mt_n = phi_i*(t_on);
534 *   double kt_n = grad_phi_i*(grad_t_on);
535 *   double mpt_n = phi_i*(latent_heat*p_on);
536 *   double mt_t = phi_i*(t_ot);
537 *   double kt_t = grad_phi_i*(grad_t_ot);
538 *   double mpt_t = phi_i*(latent_heat*p_ot);
539 *  
540 *   double f2n = mt_n + time_step*theta*kt_n - mpt_n - mt_t + time_step*(1-theta)*kt_t + mpt_t; //f2 at last newton iteration
541 *  
542 * @endcode
543 *
544 * Assembling RHS vector
545 *
546 * @code
547 *   cell_rhs(i) -= (f1n + f2n)*fe_values.JxW(q);
548 *   }
549 *   }
550 *  
551 *   cell->get_dof_indices(local_dof_indices);
552 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
553 *   {
554 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
555 *   jacobian_matrix.add(local_dof_indices[i],
556 *   local_dof_indices[j],
557 *   cell_matrix(i, j));
558 *   system_rhs(local_dof_indices[i]) += cell_rhs(i);
559 *   }
560 *   }
561 *   }
562 *  
563 *   jacobian_matrix.compress(VectorOperation::add);
564 *   system_rhs.compress(VectorOperation::add);
565 *  
566 * @endcode
567 *
568 * Applying zero BC
569 *
570 * @code
571 *   std::map<types::global_dof_index, double> boundary_values;
573 *   1,
575 *   boundary_values);
576 *   MatrixTools::apply_boundary_values(boundary_values,
577 *   jacobian_matrix,
578 *   solution_update,
579 *   system_rhs, false);
580 *  
581 *   jacobian_matrix.compress(VectorOperation::insert);
582 *  
583 *   system_rhs.compress(VectorOperation::insert);
584 *   }
585 * @endcode
586
587
588<a name="ann-get_random_number.cpp"></a>
589<h1>Annotated version of get_random_number.cpp</h1>
590 *
591 *
592 *
593 *
594 * @code
595 *   /* -----------------------------------------------------------------------------
596 *   *
597 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
598 *   * Copyright (C) 2024 by Umair Hussain
599 *   *
600 *   * This file is part of the deal.II code gallery.
601 *   *
602 *   * -----------------------------------------------------------------------------
603 *   */
604 *  
605 *   #include "PhaseFieldSolver.h"
606 *   #include <random>
607 *  
608 *   float PhaseFieldSolver::get_random_number()
609 *   {
610 *   static std::default_random_engine e;
611 *   static std::uniform_real_distribution<> dis(-0.5, 0.5); // returns a random number in the range of -0.5 to 0.5
612 *   return dis(e);
613 *   }
614 * @endcode
615
616
617<a name="ann-grid_dof.cpp"></a>
618<h1>Annotated version of grid_dof.cpp</h1>
619 *
620 *
621 *
622 *
623 * @code
624 *   /* -----------------------------------------------------------------------------
625 *   *
626 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
627 *   * Copyright (C) 2024 by Umair Hussain
628 *   *
629 *   * This file is part of the deal.II code gallery.
630 *   *
631 *   * -----------------------------------------------------------------------------
632 *   */
633 *  
634 *   #include "PhaseFieldSolver.h"
635 *  
636 *   void PhaseFieldSolver::make_grid_and_dofs() {
637 * @endcode
638 *
639 * Reading mesh
640 *
641 * @code
642 *   gridin.attach_triangulation(triangulation);
643 *   std::ifstream f("mesh/Kobayashi_mesh100x400.msh");
644 *   gridin.read_msh(f);
645 *  
646 *  
648 *   dof_handler.distribute_dofs(fe);
649 *   DoFRenumbering::subdomain_wise(dof_handler);
650 *  
651 *   DynamicSparsityPattern dsp(dof_handler.n_dofs());
652 *   DoFTools::make_sparsity_pattern(dof_handler, dsp);
653 *  
654 *   const std::vector<IndexSet> locally_owned_dofs_per_proc =
656 *   const IndexSet locally_owned_dofs =
657 *   locally_owned_dofs_per_proc[this_mpi_process];
658 *   jacobian_matrix.reinit(locally_owned_dofs,
659 *   locally_owned_dofs,
660 *   dsp,
661 *   mpi_communicator);
662 *   old_solution.reinit(locally_owned_dofs, mpi_communicator);
663 *   system_rhs.reinit(locally_owned_dofs, mpi_communicator);
664 *   conv_solution.reinit(locally_owned_dofs, mpi_communicator);
665 *   solution_update.reinit(locally_owned_dofs, mpi_communicator);
666 *  
667 *   conv_solution_np.reinit(dof_handler.n_dofs());
668 *   old_solution_np.reinit(dof_handler.n_dofs());
669 *   }
670 * @endcode
671
672
673<a name="ann-main.cpp"></a>
674<h1>Annotated version of main.cpp</h1>
675 *
676 *
677 *
678 *
679 * @code
680 *   /* -----------------------------------------------------------------------------
681 *   *
682 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
683 *   * Copyright (C) 2024 by Umair Hussain
684 *   *
685 *   * This file is part of the deal.II code gallery.
686 *   *
687 *   * -----------------------------------------------------------------------------
688 *   */
689 *  
690 *   #include <iostream>
691 *   #include "PhaseFieldSolver.h"
692 *  
693 *   int main(int argc, char **argv) {
694 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
695 *   PhaseFieldSolver phasefieldsolver;
696 *   phasefieldsolver.run();
697 *   return 0;
698 *   }
699 * @endcode
700
701
702<a name="ann-output_results.cpp"></a>
703<h1>Annotated version of output_results.cpp</h1>
704 *
705 *
706 *
707 *
708 * @code
709 *   /* -----------------------------------------------------------------------------
710 *   *
711 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
712 *   * Copyright (C) 2024 by Umair Hussain
713 *   *
714 *   * This file is part of the deal.II code gallery.
715 *   *
716 *   * -----------------------------------------------------------------------------
717 *   */
718 *  
719 *   #include "PhaseFieldSolver.h"
720 *  
721 *   void PhaseFieldSolver::output_results(const unsigned int timestep_number) const {
722 *   const Vector<double> localized_solution(old_solution);
723 *  
724 * @endcode
725 *
726 * using only one process to output the result
727 *
728 * @code
729 *   if (this_mpi_process == 0)
730 *   {
731 *   DataOut<2> data_out;
732 *   data_out.attach_dof_handler(dof_handler);
733 *  
734 *   std::vector<std::string> solution_names;
735 *   solution_names.emplace_back ("p");
736 *   solution_names.emplace_back ("T");
737 *  
738 *   data_out.add_data_vector(localized_solution, solution_names);
739 *   const std::string filename =
740 *   "solution-" + Utilities::int_to_string(timestep_number, 3) + ".vtk";
741 *   DataOutBase::VtkFlags vtk_flags;
742 *   vtk_flags.compression_level =
744 *   data_out.set_flags(vtk_flags);
745 *   std::ofstream output(filename);
746 *  
747 *   data_out.build_patches();
748 *   data_out.write_vtk(output);
749 *   }
750 *   }
751 * @endcode
752
753
754<a name="ann-run.cpp"></a>
755<h1>Annotated version of run.cpp</h1>
756 *
757 *
758 *
759 *
760 * @code
761 *   /* -----------------------------------------------------------------------------
762 *   *
763 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
764 *   * Copyright (C) 2024 by Umair Hussain
765 *   *
766 *   * This file is part of the deal.II code gallery.
767 *   *
768 *   * -----------------------------------------------------------------------------
769 *   */
770 *  
771 *   #include "PhaseFieldSolver.h"
772 *   #include <time.h>
773 *  
774 *   void PhaseFieldSolver::run() {
775 *   make_grid_and_dofs();
776 *   pcout << "Processors used: " << n_mpi_processes << std::endl;
777 *   pcout << " Number of degrees of freedom: " << dof_handler.n_dofs()
778 *   << " (by partition:";
779 *   for (unsigned int p = 0; p < n_mpi_processes; ++p)
780 *   pcout << (p == 0 ? ' ' : '+')
782 *   p));
783 *   pcout << ")" << std::endl;
784 * @endcode
785 *
786 * Initialise the solution
787 *
788 * @code
789 *   InitialValues initial_value;
790 *   VectorTools::interpolate(dof_handler,
791 *   initial_value,
792 *   old_solution);
793 *   VectorTools::interpolate(dof_handler,
794 *   initial_value,
795 *   conv_solution);
796 * @endcode
797 *
798 * Applying Boundary Conditions at t=0
799 *
800 * @code
801 *   applying_bc();
802 * @endcode
803 *
804 * Plotting initial solution
805 *
806 * @code
807 *   output_results(0);
808 *  
809 * @endcode
810 *
811 * Time steps begin here:
812 *
813 * @code
814 *   unsigned int timestep_number = 1;
815 *   for (; time <= final_time; time += time_step, ++timestep_number) {
816 *  
817 *   pcout << "Time step " << timestep_number << " at t=" << time+time_step
818 *   << std::endl;
819 *  
820 *   conv_solution.operator=(old_solution); // initialising the newton solution
821 *  
822 * @endcode
823 *
824 * Newton-Raphson iterations begin here:
825 *
826 * @code
827 *   for (unsigned int it = 1; it <= 100; ++it) {
828 *   pcout << "Newton iteration number:" << it << std::endl;
829 *  
830 *   if (it == 100) {
831 *   pcout << "Convergence Failure!!!!!!!!!!!!!!!" << std::endl;
832 *   std::exit(0);
833 *   }
834 * @endcode
835 *
836 * Saving parallel vectors as non-parallel ones
837 *
838 * @code
839 *   conv_solution_np = conv_solution;
840 *   old_solution_np = old_solution;
841 * @endcode
842 *
843 * Initialise the delta solution as zero
844 *
845 * @code
846 *   VectorTools::interpolate(dof_handler,
848 *   solution_update);
849 *   solution_update.compress(VectorOperation::insert);
850 * @endcode
851 *
852 * Assemble Jacobian and Residual
853 *
854 * @code
855 *   assemble_system();
856 * @endcode
857 *
858 * Solving to get delta solution
859 *
860 * @code
861 *   solve();
862 * @endcode
863 *
864 * Checking for convergence
865 *
866 * @code
867 *   double residual_norm = system_rhs.l2_norm(); //the norm of residual should converge to zero as the solution converges
868 * @endcode
869 *
870 * pcout << "Nothing wrong till here!!!!!!" << std::endl;
871 *
872 * @code
873 *   pcout << "the residual is:" << residual_norm << std::endl;
874 *   if (residual_norm <= (1e-4)) {
875 *   pcout << "Solution Converged!" << std::endl;
876 *   break; //Break to next time step if the N-R iterations converge
877 *   }
878 *   }
879 * @endcode
880 *
881 * Transfer the converged solution to the old_solution vector to plot output
882 *
883 * @code
884 *   old_solution.operator=(conv_solution);
885 *   old_solution.compress(VectorOperation::insert);
886 * @endcode
887 *
888 * output the solution at only specific number of time steps
889 *
890 * @code
891 *   if (timestep_number%10 == 0)
892 *   output_results(timestep_number);
893 *   }
894 *   }
895 * @endcode
896
897
898<a name="ann-solve.cpp"></a>
899<h1>Annotated version of solve.cpp</h1>
900 *
901 *
902 *
903 *
904 * @code
905 *   /* -----------------------------------------------------------------------------
906 *   *
907 *   * SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
908 *   * Copyright (C) 2024 by Umair Hussain
909 *   *
910 *   * This file is part of the deal.II code gallery.
911 *   *
912 *   * -----------------------------------------------------------------------------
913 *   */
914 *  
915 *   #include "PhaseFieldSolver.h"
916 *  
917 *   void PhaseFieldSolver::solve(){
918 * @endcode
919 *
920 * Using a direct parallel solver
921 *
922 * @code
923 *   SolverControl cn;
924 *   PETScWrappers::SparseDirectMUMPS A_direct(cn, mpi_communicator);
925 *   A_direct.solve(jacobian_matrix, solution_update, system_rhs);
926 * @endcode
927 *
928 * Updating the solution by adding the delta solution
929 *
930 * @code
931 *   conv_solution.add(1, solution_update);
932 *   conv_solution.compress(VectorOperation::add);
933 *   }
934 * @endcode
935
936
937*/
std::vector< bool > component_mask
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
Definition fe_q.h:554
Definition point.h:111
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
const Event initial
Definition event.cc:64
void subdomain_wise(DoFHandler< dim, spacedim > &dof_handler)
void random(DoFHandler< dim, spacedim > &dof_handler)
std::vector< IndexSet > locally_owned_dofs_per_subdomain(const DoFHandler< dim, spacedim > &dof_handler)
unsigned int count_dofs_with_subdomain_association(const DoFHandler< dim, spacedim > &dof_handler, const types::subdomain_id subdomain)
void partition_triangulation(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
@ matrix
Contents is actually a matrix.
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:74
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
void apply(const Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space >::member_type &team_member, const Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_data, const ViewTypeIn in, ViewTypeOut out)
VectorType::value_type * begin(VectorType &V)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:92
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:107
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
void interpolate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={}, const unsigned int level=numbers::invalid_unsigned_int)
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask={})
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
long double gamma(const unsigned int n)
void copy(const T *begin, const T *end, U *dest)
int(&) functions(const void *v1, const void *v2)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
inline ::VectorizedArray< Number, width > atan(const ::VectorizedArray< Number, width > &x)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
DataOutBase::CompressionLevel compression_level