void Diffusion::flux_matrix(const FEMesh *mesh, const Element *el, const ElementFuncNodeIterator &j, const Flux *flux, const MasterPosition &pt, double time, SmallSystem *fluxdata) const { // The atom flux matrix M_{ij} multiplies the vector of nodal // concentrations to give the vector atom current J at point pt. // M_{ij} = -D_{ik} grad_k N_j // J_i = -M_{ij} T_j // where N_j is the shapefunction at node j. // printf("flux"); if (*flux != *atom_flux) { throw ErrProgrammingError("Unexpected flux", __FILE__, __LINE__); } double sf = j.shapefunction( pt ); double dsf0 = j.dshapefunction( 0, pt ); double dsf1 = j.dshapefunction( 1, pt ); #if DIM==3 double dsf2 = j.dshapefunction( 2, pt ); #endif const SymmMatrix3 cond( conductivitytensor( mesh, el, pt ) ); // Loop over flux components. Loop over all components, even if // the flux is in-plane, because the out-of-plane components of // the flux matrix are used to construct the constraint equation. for(VectorFieldIterator i; !i.end(); ++i){ #if DIM==2 // in-plane concentration gradient contributions fluxdata->stiffness_matrix_element( i, concentration, j ) -= cond(i.integer(), 0) * dsf0 + cond(i.integer(), 1) * dsf1; // out-of-plane concentration gradient contribution if(!concentration->in_plane(mesh)) fluxdata->stiffness_matrix_element(i, concentration->out_of_plane(), j) -= cond(i.integer(), 2) * sf; #elif DIM==3 fluxdata->stiffness_matrix_element( i, concentration, j ) -= cond( i.integer(), 0 ) * dsf0 + cond( i.integer(), 1 ) * dsf1 + cond( i.integer(), 2 ) * dsf2; #endif } } // end of 'Diffusion::flux_matrix'
void CLargeStrainElasticity::flux_matrix(const FEMesh *mesh, const Element *element, const ElementFuncNodeIterator &node, const Flux *flux, const MasterPosition &pt, double time, SmallSystem *fluxmtx) const { int (*ij2voigt)(int,int) = &SymTensorIndex::ij2voigt; // shorter func name SmallMatrix dU(3); // gradient of displacement double Fval; // the value of the shape function (for node) DoubleVec dF(3); // and its derivative at the given pt bool inplane = false; // needed both in 2D & 3D versions regardless, // passed to contract_C_dU_dF #if DIM==2 // in 2D, check if it is an in-plane eqn or a plane-flux eqn. static CompoundField *displacement = dynamic_cast<CompoundField*>(Field::getField("Displacement")); inplane = displacement->in_plane( mesh ); #endif // check for unexpected flux, flux should be a stress flux if (*flux != *stress_flux) { throw ErrProgrammingError("Unexpected flux", __FILE__, __LINE__); } // evaluate the shape function and its gradient (of node) at the given pt Fval = node.shapefunction( pt ); // value of the shape function dF[0] = node.dshapefunction( 0, pt ); // x-deriv of the shape function dF[1] = node.dshapefunction( 1, pt ); // y-deriv of the shape function #if DIM==3 dF[2] = node.dshapefunction( 2, pt ); // z-deriv of the shape function #endif computeDisplacementGradient( mesh, element, pt, dU ); const Cijkl CC = cijkl( mesh, element, pt ); // elasticity modulus // add the flux contributions to stiffness matrix element // k_indx is needed for fluxmtx->stifness_matrix_element function, // which does not take int k as argument VectorFieldIndex k_indx; for (SymTensorIterator ij_iter; !ij_iter.end(); ++ij_iter) { int k0, k1, k2, ij = ij_iter.integer(); double nonlinear_part; // to store the sum from the nonlinear terms // TODO: Use tensor iterators for k0, k1, k2. #if DIM==2 // sum CC(i,j,k,l)*dF(l), k=0 over l=0,1, then add to stiffness_mtx k_indx.set( 0 ); k0 = ij2voigt( 0,0 ); k1 = ij2voigt( 0,1 ); nonlinear_part = contract_C_dU_dF(CC, dU, dF, ij, 0, inplane ); // at ij, k=0 fluxmtx->stiffness_matrix_element( ij_iter, displacement, k_indx, node ) += CC( ij,k0 ) * dF[0] + CC( ij,k1 ) * dF[1] + nonlinear_part; // sum CC(i,j,k,l)*dF(l), k=1 over l=0,1, then add to stiffness_mtx k_indx.set( 1 ); k0 = ij2voigt( 1,0 ); k1 = ij2voigt( 1,1 ); nonlinear_part = contract_C_dU_dF( CC, dU, dF, ij, 1, inplane ); // at ij, k=1 fluxmtx->stiffness_matrix_element( ij_iter, displacement, k_indx, node ) += CC( ij,k0 ) * dF[0] + CC( ij,k1 ) * dF[1] + nonlinear_part; #elif DIM==3 // sum CC(i,j,k,l)*dF(l), k=0 over l=0,1,2 then add to stiffness_mtx k_indx.set( 0 ); k0 = ij2voigt( 0,0 ); k1 = ij2voigt( 0,1 ); k2 = ij2voigt( 0,2 ); nonlinear_part = contract_C_dU_dF( CC, dU, dF, ij, 0, inplane ); // at ij, k=0 fluxmtx->stiffness_matrix_element( ij_iter, displacement, k_indx, node ) += CC( ij,k0 ) * dF[0] + CC( ij,k1 ) * dF[1] + CC( ij,k2 ) * dF[2] + nonlinear_part; // sum CC(i,j,k,l)*dF(l), k=1 over l=0,1,2 then add to stiffness_mtx k_indx.set( 1 ); k0 = ij2voigt( 1,0 ); k1 = ij2voigt( 1,1 ); k2 = ij2voigt( 1,2 ); nonlinear_part = contract_C_dU_dF( CC, dU, dF, ij, 1, inplane ); // at ij, k=1 fluxmtx->stiffness_matrix_element( ij_iter, displacement, k_indx, node ) += CC( ij,k0 ) * dF[0] + CC( ij,k1 ) * dF[1] + CC( ij,k2 ) * dF[2] + nonlinear_part; // sum CC(i,j,k,l)*dF(l), k=2 over l=0,1,2 then add to stiffness_mtx k_indx.set( 2 ); k0 = ij2voigt( 2,0 ); k1 = ij2voigt( 2,1 ); k2 = ij2voigt( 2,2 ); nonlinear_part = contract_C_dU_dF( CC, dU, dF, ij, 2, inplane ); // at ij, k=2 fluxmtx->stiffness_matrix_element( ij_iter, displacement, k_indx, node ) += CC( ij,k0 ) * dF[0] + CC( ij,k1 ) * dF[1] + CC( ij,k2 ) * dF[2] + nonlinear_part; #endif #if DIM==2 if ( !inplane ) // now contributions from z-deriv of displacement field { Field *disp_z_deriv = displacement->out_of_plane(); for(IteratorP k_iter = disp_z_deriv->iterator( ALL_INDICES ); !k_iter.end(); ++k_iter) { double diag_factor = ( k_iter.integer()==2 ? 1.0 : 0.5 ); k2 = ij2voigt( 2, k_iter.integer() ); fluxmtx->stiffness_matrix_element( ij_iter, disp_z_deriv, k_iter, node ) += diag_factor * Fval * CC( ij,k2 ); } } // end of 'if (!inplane)' #endif } // end of loop over ij } // end of 'CLargeStrainElasticity::flux_matrix'