v0.14.0
initial_diffusion.cpp
Go to the documentation of this file.
1 /**
2  * \file inital_diffusion.cpp
3  * \example inital_diffusion.cpp
4  *
5  **/
6 
7 #include <stdlib.h>
8 #include <cmath>
9 #include <MoFEM.hpp>
10 
11 #define BOOST_MATH_GAUSS_NO_COMPUTE_ON_DEMAND
12 
13 using namespace MoFEM;
14 
15 static char help[] = "...\n\n";
16 
17 template <int DIM> struct ElementsAndOps {};
18 
19 //! [Define dimension]
20 constexpr int SPACE_DIM = 3; //< Space dimension of problem, mesh
21 //! [Define dimension]
22 
26 
32  PETSC>::LinearForm<GAUSS>::OpBaseTimesScalar<1>;
36  PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
37 
38 const double n = 1.44; ///< refractive index of diffusive medium
39 const double c = 30.; ///< speed of light (cm/ns)
40 const double v = c / n; ///< phase velocity of light in medium (cm/ns)
41 
42 double mu_a; ///< absorption coefficient (cm^-1)
43 double mu_sp; ///< scattering coefficient (cm^-1)
44 double D;
45 
47 double beam_radius; //< spot radius
50 double flux_magnitude = 1e3; ///< impulse magnitude
51 double initial_time;
52 
53 char out_file_name[255] = "init_file.dat";;
54 int numHoLevels = 1;
55 
56 PetscBool output_volume = PETSC_FALSE;
57 
58 #include <boost/math/quadrature/gauss_kronrod.hpp>
59 using namespace boost::math::quadrature;
60 
61 struct OpError : public DomainEleOp {
62  OpError(boost::shared_ptr<VectorDouble> u_at_pts_ptr, double &l2_error)
63  : DomainEleOp("PHOTON_FLUENCE_RATE", OPROW), uAtPtsPtr(u_at_pts_ptr),
64  l2Error(l2_error) {
65  // Only will be executed once on vertices
66  std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
67  }
68 
69 private:
70  boost::shared_ptr<VectorDouble> uAtPtsPtr;
71  double &l2Error;
72 };
73 
75 public:
77 
78  // Declaration of the main function to run analysis
79  MoFEMErrorCode runProgram();
80 
81  /**
82  * @brief Pulse is infinitely short.
83  *
84  * \note It is good approximation of pulse in femtosecond scale. To make it
85  * longer one can apply third integral over time.
86  *
87  * \note Note analysis in photon_diffusion is shifted in time by initial_time.
88  *
89  * @param x
90  * @param y
91  * @param z
92  * @return double
93  */
94  static double sourceFunction(const double x, const double y, const double z) {
95  const double A = 4. * D * v * initial_time;
96  const double T =
97  (v / pow(M_PI * A, 3. / 2.)) * exp(-mu_a * v * initial_time);
98 
99  auto phi_pulse = [&](const double r_s, const double phi_s) {
100  const double xs = r_s * cos(phi_s);
101  const double ys = r_s * sin(phi_s);
102  const double xp = x - xs - beam_centre_x;
103  const double yp = y - ys - beam_centre_y;
104  const double zp1 = z + slab_thickness / 2. - 1. / mu_sp;
105  const double zp2 = z + slab_thickness / 2. + 1. / mu_sp;
106  const double P1 = xp * xp + yp * yp + zp1 * zp1;
107  const double P2 = xp * xp + yp * yp + zp2 * zp2;
108  return r_s * (exp(-P1 / A) - exp(-P2 / A));
109  };
110 
111  auto f = [&](const double r_s) {
112  auto g = [&](const double phi_s) { return phi_pulse(r_s, phi_s); };
113  return gauss_kronrod<double, 15>::integrate(
114  g, 0, 2 * M_PI, 0, std::numeric_limits<float>::epsilon());
115  };
116 
117  return T * flux_magnitude *
118  gauss_kronrod<double, 15>::integrate(
119  f, 0, beam_radius, 0, std::numeric_limits<float>::epsilon());
120  };
121 
122 private:
123  // Declaration of other main functions called in runProgram()
124  MoFEMErrorCode readMesh();
125  MoFEMErrorCode setupProblem();
126  MoFEMErrorCode setIntegrationRules();
127  MoFEMErrorCode initialCondition();
128  MoFEMErrorCode boundaryCondition();
129  MoFEMErrorCode assembleSystem();
130  MoFEMErrorCode solveSystem();
131  MoFEMErrorCode checkResults();
132  MoFEMErrorCode outputResults();
133 
134  // Main interfaces
136 };
137 
139 
142 
143  auto *simple = mField.getInterface<Simple>();
145  CHKERR simple->getOptions();
146  CHKERR simple->loadFile();
147 
149 }
150 
153 
154  auto *simple = mField.getInterface<Simple>();
155  CHKERR simple->addDomainField("PHOTON_FLUENCE_RATE", H1,
157  CHKERR simple->addBoundaryField("PHOTON_FLUENCE_RATE", H1,
159 
160  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-flux_magnitude",
161  &flux_magnitude, PETSC_NULL);
162  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-slab_thickness",
163  &slab_thickness, PETSC_NULL);
164  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-beam_radius", &beam_radius,
165  PETSC_NULL);
166  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-beam_centre_x", &beam_centre_x,
167  PETSC_NULL);
168  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-beam_centre_y", &beam_centre_y,
169  PETSC_NULL);
170  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-mu_a", &mu_a, PETSC_NULL);
171  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-mu_sp", &mu_sp, PETSC_NULL);
172  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-initial_time", &initial_time,
173  PETSC_NULL);
174 
175  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-output_file", out_file_name,
176  255, PETSC_NULL);
177  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-output_volume", &output_volume,
178  PETSC_NULL);
179 
180  D = 1. / (3. * (mu_a + mu_sp));
181 
182  MOFEM_LOG("INITIAL", Sev::inform) << "Refractive index: " << n;
183  MOFEM_LOG("INITIAL", Sev::inform) << "Speed of light (cm/ns): " << c;
184  MOFEM_LOG("INITIAL", Sev::inform)
185  << "Phase velocity in medium (cm/ns): " << v;
186  MOFEM_LOG("INITIAL", Sev::inform)
187  << "Absorption coefficient (cm^-1): " << mu_a;
188  MOFEM_LOG("INITIAL", Sev::inform)
189  << "Scattering coefficient (cm^-1): " << mu_sp;
190  MOFEM_LOG("INITIAL", Sev::inform) << "Diffusion coefficient D : " << D;
191  MOFEM_LOG("INITIAL", Sev::inform) << "Impulse magnitude: " << flux_magnitude;
192  MOFEM_LOG("INITIAL", Sev::inform) << "Compute time (ns): " << initial_time;
193  MOFEM_LOG("INITIAL", Sev::inform) << "Slab thickness: " << slab_thickness;
194 
195  int order = 2;
196  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
197 
198  MOFEM_LOG("INITIAL", Sev::inform) << "Approximation order: " << order;
199 
200  CHKERR simple->setFieldOrder("PHOTON_FLUENCE_RATE", order);
201 
202  // if (numHoLevels > 0) {
203 
204  // Range ho_ents;
205  // for (_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(mField, BLOCKSET, it)) {
206  // if (it->getName().compare(0, 3, "CAM") == 0) {
207  // CHKERR mField.get_moab().get_entities_by_dimension(it->getMeshset(), 2,
208  // ho_ents, true);
209  // }
210  // }
211 
212  // EntityHandle meshset;
213  // CHKERR mField.get_moab().create_meshset(MESHSET_SET, meshset);
214  // CHKERR mField.get_moab().add_entities(meshset, ho_ents);
215  // std::string field_name;
216  // field_name = "out_test_" +
217  // boost::lexical_cast<std::string>(mField.get_comm_rank()) +
218  // ".vtk";
219  // CHKERR mField.get_moab().write_file(field_name.c_str(), "VTK", "", &meshset,
220  // 1);
221  // CHKERR mField.get_moab().delete_entities(&meshset, 1);
222 
223  // CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(ho_ents);
224 
225  // CHKERR simple->setFieldOrder("PHOTON_FLUENCE_RATE", order + 1, &ho_ents);
226 
227  // CHKERR mField.getInterface<CommInterface>()->synchroniseFieldEntities(
228  // "PHOTON_FLUENCE_RATE");
229  // }
230 
231  CHKERR simple->setUp();
232 
234 }
235 
238 
239  auto integration_rule = [](int o_row, int o_col, int approx_order) {
240  return 2 * approx_order;
241  };
242 
243  auto *pipeline_mng = mField.getInterface<PipelineManager>();
244  CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
245  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
246  CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
247  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
248 
250 }
251 
255 }
256 
259 
260  auto *simple = mField.getInterface<Simple>();
261 
262  // Get boundary edges marked in block named "BOUNDARY_CONDITION"
263  Range boundary_ents;
265  std::string entity_name = it->getName();
266  if (entity_name.compare(0, 3, "INT") == 0) {
267  CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 1,
268  boundary_ents, true);
269  }
270  }
271  // Add vertices to boundary entities
272  Range boundary_verts;
273  CHKERR mField.get_moab().get_connectivity(boundary_ents, boundary_verts,
274  true);
275  boundary_ents.merge(boundary_verts);
276 
277  // Remove DOFs as homogeneous boundary condition is used
278  CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
279  simple->getProblemName(), "PHOTON_FLUENCE_RATE", boundary_ents);
280 
282 }
283 
286 
287  auto bc_mng = mField.getInterface<BcManager>();
288  auto &bc_map = bc_mng->getBcMapByBlockName();
289 
290  auto add_domain_base_ops = [&](auto &pipeline) {
291  auto jac_ptr = boost::make_shared<MatrixDouble>();
292  auto det_ptr = boost::make_shared<VectorDouble>();
293  pipeline.push_back(new OpCalculateHOJacVolume(jac_ptr));
294  pipeline.push_back(new OpInvertMatrix<3>(jac_ptr, det_ptr, nullptr));
295  pipeline.push_back(new OpSetHOWeights(det_ptr));
296  };
297 
298  auto add_domain_lhs_ops = [&](auto &pipeline) {
299  pipeline.push_back(new OpDomainMass(
300  "PHOTON_FLUENCE_RATE", "PHOTON_FLUENCE_RATE",
301  [](const double, const double, const double) { return 1; }));
302  };
303 
304  auto add_domain_rhs_ops = [&](auto &pipeline) {
305  pipeline.push_back(
306  new OpDomainSource("PHOTON_FLUENCE_RATE", sourceFunction));
307  };
308 
309  auto pipeline_mng = mField.getInterface<PipelineManager>();
310  add_domain_base_ops(pipeline_mng->getOpDomainLhsPipeline());
311  add_domain_base_ops(pipeline_mng->getOpDomainRhsPipeline());
312  add_domain_lhs_ops(pipeline_mng->getOpDomainLhsPipeline());
313  add_domain_rhs_ops(pipeline_mng->getOpDomainRhsPipeline());
314 
316 }
317 
320  auto *simple = mField.getInterface<Simple>();
321  auto *pipeline_mng = mField.getInterface<PipelineManager>();
322  auto solver = pipeline_mng->createKSP();
323 
324  CHKERR KSPSetFromOptions(solver);
325  // CHKERR KSPSetUp(solver);
326 
327  auto dm = simple->getDM();
328  auto X = createDMVector(dm);
329  auto F = vectorDuplicate(X);
330 
331  MOFEM_LOG("INITIAL", Sev::inform) << "Solver start";
332  CHKERR KSPSolve(solver, F, X);
333  CHKERR VecGhostUpdateBegin(X, INSERT_VALUES, SCATTER_FORWARD);
334  CHKERR VecGhostUpdateEnd(X, INSERT_VALUES, SCATTER_FORWARD);
335  CHKERR DMoFEMMeshToLocalVector(dm, X, INSERT_VALUES, SCATTER_REVERSE);
336 
337  MOFEM_LOG("INITIAL", Sev::inform)
338  << "writing vector in binary to " << out_file_name << " ...";
339  PetscViewer viewer;
340  PetscViewerBinaryOpen(PETSC_COMM_WORLD, out_file_name, FILE_MODE_WRITE,
341  &viewer);
342  VecView(X, viewer);
343  PetscViewerDestroy(&viewer);
344 
345  MOFEM_LOG("INITIAL", Sev::inform) << "Solver done";
347 }
348 
351  auto *pipeline_mng = mField.getInterface<PipelineManager>();
352  pipeline_mng->getDomainLhsFE().reset();
353  auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
354 
355  auto u_ptr = boost::make_shared<VectorDouble>();
356  post_proc_fe->getOpPtrVector().push_back(
357  new OpCalculateScalarFieldValues("U", u_ptr));
358 
360 
361  post_proc_fe->getOpPtrVector().push_back(
362 
363  new OpPPMap(
364 
365  post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
366 
367  {{"PHOTON_FLUENCE_RATE", u_ptr}},
368 
369  {},
370 
371  {},
372 
373  {})
374 
375  );
376 
377  pipeline_mng->getDomainRhsFE() = post_proc_fe;
378  CHKERR pipeline_mng->loopFiniteElements();
379  CHKERR post_proc_fe->writeFile("out_initial.h5m");
381 }
382 
385 
386  CHKERR readMesh();
393  if (output_volume)
395 
397 }
398 
399 int main(int argc, char *argv[]) {
400 
401  // Initialisation of MoFEM/PETSc and MOAB data structures
402  const char param_file[] = "param_file.petsc";
403  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
404 
405  // Add logging channel for example
406  auto core_log = logging::core::get();
407  core_log->add_sink(
408  LogManager::createSink(LogManager::getStrmWorld(), "INITIAL"));
409  LogManager::setLog("INITIAL");
410  MOFEM_LOG_TAG("INITIAL", "initial_diffusion")
411 
412  // Error handling
413  try {
414  // Register MoFEM discrete manager in PETSc
415  DMType dm_name = "DMMOFEM";
416  CHKERR DMRegister_MoFEM(dm_name);
417 
418  // Create MOAB instance
419  moab::Core mb_instance; // mesh database
420  moab::Interface &moab = mb_instance; // mesh database interface
421 
422  // Create MoFEM instance
423  MoFEM::Core core(moab); // finite element database
424  MoFEM::Interface &m_field = core; // finite element interface
425 
426  // Run the main analysis
427  PhotonDiffusion heat_problem(m_field);
428  CHKERR heat_problem.runProgram();
429  }
430  CATCH_ERRORS;
431 
432  // Finish work: cleaning memory, getting statistics, etc.
434 
435  return 0;
436 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
g
constexpr double g
Definition: shallow_wave.cpp:63
beam_centre_x
double beam_centre_x
Definition: initial_diffusion.cpp:48
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
OpError::l2Error
double & l2Error
Definition: initial_diffusion.cpp:71
OpError
Definition: initial_diffusion.cpp:61
mu_sp
double mu_sp
scattering coefficient (cm^-1)
Definition: initial_diffusion.cpp:43
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, SPACE_DIM > OpDomainGradTimesVec
Definition: initial_diffusion.cpp:34
OpError::uAtPtsPtr
boost::shared_ptr< VectorDouble > uAtPtsPtr
Definition: initial_diffusion.cpp:70
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
PhotonDiffusion
Definition: initial_diffusion.cpp:74
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
beam_radius
double beam_radius
Definition: initial_diffusion.cpp:47
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpDomainMass
Definition: initial_diffusion.cpp:28
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
PhotonDiffusion::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
Definition: initial_diffusion.cpp:236
beam_centre_y
double beam_centre_y
Definition: initial_diffusion.cpp:49
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:523
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
out_file_name
char out_file_name[255]
Definition: initial_diffusion.cpp:53
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
flux_magnitude
double flux_magnitude
impulse magnitude
Definition: initial_diffusion.cpp:50
PhotonDiffusion::PhotonDiffusion
PhotonDiffusion(MoFEM::Interface &m_field)
Definition: initial_diffusion.cpp:138
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Definition: child_and_parent.cpp:53
PhotonDiffusion::readMesh
MoFEMErrorCode readMesh()
Definition: initial_diffusion.cpp:140
order
constexpr int order
Definition: dg_projection.cpp:18
slab_thickness
double slab_thickness
Definition: initial_diffusion.cpp:46
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
PhotonDiffusion::setupProblem
MoFEMErrorCode setupProblem()
Definition: initial_diffusion.cpp:151
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
OpDomainTimesScalarField
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
Definition: initial_diffusion.cpp:32
D
double D
Definition: initial_diffusion.cpp:44
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
PhotonDiffusion::mField
MoFEM::Interface & mField
Definition: initial_diffusion.cpp:135
PhotonDiffusion::initialCondition
MoFEMErrorCode initialCondition()
Definition: initial_diffusion.cpp:252
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
mu_a
double mu_a
absorption coefficient (cm^-1)
Definition: initial_diffusion.cpp:42
MoFEM::BcManager::getBcMapByBlockName
BcMapByBlockName & getBcMapByBlockName()
Get the bc map.
Definition: BcManager.hpp:243
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
PhotonDiffusion::runProgram
MoFEMErrorCode runProgram()
Definition: initial_diffusion.cpp:383
n
const double n
refractive index of diffusive medium
Definition: initial_diffusion.cpp:38
OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpDomainSource
Definition: initial_diffusion.cpp:36
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
PhotonDiffusion::solveSystem
MoFEMErrorCode solveSystem()
Definition: initial_diffusion.cpp:318
OpError::OpError
OpError(boost::shared_ptr< VectorDouble > u_at_pts_ptr, double &l2_error)
Definition: initial_diffusion.cpp:62
BiLinearForm
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
PhotonDiffusion::boundaryCondition
MoFEMErrorCode boundaryCondition()
Definition: initial_diffusion.cpp:257
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
help
static char help[]
Definition: initial_diffusion.cpp:15
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
ElementsAndOps
Definition: child_and_parent.cpp:18
Range
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
DomainEleOp
PhotonDiffusion::assembleSystem
MoFEMErrorCode assembleSystem()
Definition: initial_diffusion.cpp:284
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: helmholtz.cpp:25
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
initial_time
double initial_time
Definition: initial_diffusion.cpp:51
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
approx_order
int approx_order
Definition: test_broken_space.cpp:50
output_volume
PetscBool output_volume
Definition: initial_diffusion.cpp:56
MoFEM::OpSetHOWeights
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:159
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::OpInvertMatrix
Definition: UserDataOperators.hpp:3249
numHoLevels
int numHoLevels
Definition: initial_diffusion.cpp:54
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
Definition: child_and_parent.cpp:55
PhotonDiffusion::sourceFunction
static double sourceFunction(const double x, const double y, const double z)
Pulse is infinitely short.
Definition: initial_diffusion.cpp:94
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
SPACE_DIM
constexpr int SPACE_DIM
[Define dimension]
Definition: initial_diffusion.cpp:20
main
int main(int argc, char *argv[])
Definition: initial_diffusion.cpp:399
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
F
@ F
Definition: free_surface.cpp:394
MoFEM::OpCalculateHOJacVolume
OpCalculateHOJacForVolume OpCalculateHOJacVolume
Definition: HODataOperators.hpp:34
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
PhotonDiffusion::outputResults
MoFEMErrorCode outputResults()
Definition: initial_diffusion.cpp:349